1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 //-------------------------------------------------------
18 // Implementation of the TPC tracker
20 // Origin: Marian Ivanov Marian.Ivanov@cern.ch
22 // AliTPC parallel tracker
23 //-------------------------------------------------------
28 #include "Riostream.h"
29 #include <TClonesArray.h>
31 #include <TObjArray.h>
34 #include "AliComplexCluster.h"
36 #include "AliESDkink.h"
38 #include "AliRunLoader.h"
39 #include "AliTPCClustersRow.h"
40 #include "AliTPCParam.h"
41 #include "AliTPCReconstructor.h"
42 #include "AliTPCclusterMI.h"
43 #include "AliTPCpolyTrack.h"
44 #include "AliTPCreco.h"
45 #include "AliTPCseed.h"
46 #include "AliTPCtrackerMI.h"
47 #include "TStopwatch.h"
48 #include "AliTPCReconstructor.h"
49 #include "AliESDkink.h"
51 #include "TTreeStream.h"
52 #include "AliAlignObj.h"
53 #include "AliTrackPointArray.h"
57 ClassImp(AliTPCtrackerMI)
60 class AliTPCFastMath {
63 static Double_t FastAsin(Double_t x);
65 static Double_t fgFastAsin[20000]; //lookup table for fast asin computation
68 Double_t AliTPCFastMath::fgFastAsin[20000];
69 AliTPCFastMath gAliTPCFastMath; // needed to fill the LUT
71 AliTPCFastMath::AliTPCFastMath(){
73 // initialized lookup table;
74 for (Int_t i=0;i<10000;i++){
75 fgFastAsin[2*i] = TMath::ASin(i/10000.);
76 fgFastAsin[2*i+1] = (TMath::ASin((i+1)/10000.)-fgFastAsin[2*i]);
80 Double_t AliTPCFastMath::FastAsin(Double_t x){
82 // return asin using lookup table
84 Int_t index = int(x*10000);
85 return fgFastAsin[2*index]+(x*10000.-index)*fgFastAsin[2*index+1];
88 Int_t index = int(x*10000);
89 return -(fgFastAsin[2*index]+(x*10000.-index)*fgFastAsin[2*index+1]);
95 Int_t AliTPCtrackerMI::UpdateTrack(AliTPCseed * track, Int_t accept){
97 //update track information using current cluster - track->fCurrentCluster
100 AliTPCclusterMI* c =track->fCurrentCluster;
101 if (accept>0) track->fCurrentClusterIndex1 |=0x8000; //sign not accepted clusters
103 UInt_t i = track->fCurrentClusterIndex1;
105 Int_t sec=(i&0xff000000)>>24;
106 //Int_t row = (i&0x00ff0000)>>16;
107 track->fRow=(i&0x00ff0000)>>16;
108 track->fSector = sec;
109 // Int_t index = i&0xFFFF;
110 if (sec>=fParam->GetNInnerSector()) track->fRow += fParam->GetNRowLow();
111 track->SetClusterIndex2(track->fRow, i);
112 //track->fFirstPoint = row;
113 //if ( track->fLastPoint<row) track->fLastPoint =row;
114 // if (track->fRow<0 || track->fRow>160) {
115 // printf("problem\n");
117 if (track->fFirstPoint>track->fRow)
118 track->fFirstPoint = track->fRow;
119 if (track->fLastPoint<track->fRow)
120 track->fLastPoint = track->fRow;
123 track->fClusterPointer[track->fRow] = c;
126 Double_t angle2 = track->GetSnp()*track->GetSnp();
128 //SET NEW Track Point
130 if (angle2<1) //PH sometimes angle2 is very big. To be investigated...
132 angle2 = TMath::Sqrt(angle2/(1-angle2));
133 AliTPCTrackerPoint &point =*(track->GetTrackPoint(track->fRow));
135 point.SetSigmaY(c->GetSigmaY2()/track->fCurrentSigmaY2);
136 point.SetSigmaZ(c->GetSigmaZ2()/track->fCurrentSigmaZ2);
137 point.SetErrY(sqrt(track->fErrorY2));
138 point.SetErrZ(sqrt(track->fErrorZ2));
140 point.SetX(track->GetX());
141 point.SetY(track->GetY());
142 point.SetZ(track->GetZ());
143 point.SetAngleY(angle2);
144 point.SetAngleZ(track->GetTgl());
145 if (point.fIsShared){
146 track->fErrorY2 *= 4;
147 track->fErrorZ2 *= 4;
151 Double_t chi2 = track->GetPredictedChi2(track->fCurrentCluster);
153 track->fErrorY2 *= 1.3;
154 track->fErrorY2 += 0.01;
155 track->fErrorZ2 *= 1.3;
156 track->fErrorZ2 += 0.005;
158 if (accept>0) return 0;
159 if (track->GetNumberOfClusters()%20==0){
160 // if (track->fHelixIn){
161 // TClonesArray & larr = *(track->fHelixIn);
162 // Int_t ihelix = larr.GetEntriesFast();
163 // new(larr[ihelix]) AliHelix(*track) ;
166 track->fNoCluster =0;
167 return track->Update(c,chi2,i);
172 Int_t AliTPCtrackerMI::AcceptCluster(AliTPCseed * seed, AliTPCclusterMI * cluster, Float_t factor,
173 Float_t cory, Float_t corz)
176 // decide according desired precision to accept given
177 // cluster for tracking
178 Double_t sy2=ErrY2(seed,cluster)*cory;
179 Double_t sz2=ErrZ2(seed,cluster)*corz;
180 //sy2=ErrY2(seed,cluster)*cory;
181 //sz2=ErrZ2(seed,cluster)*cory;
183 Double_t sdistancey2 = sy2+seed->GetSigmaY2();
184 Double_t sdistancez2 = sz2+seed->GetSigmaZ2();
186 Double_t rdistancey2 = (seed->fCurrentCluster->GetY()-seed->GetY())*
187 (seed->fCurrentCluster->GetY()-seed->GetY())/sdistancey2;
188 Double_t rdistancez2 = (seed->fCurrentCluster->GetZ()-seed->GetZ())*
189 (seed->fCurrentCluster->GetZ()-seed->GetZ())/sdistancez2;
191 Double_t rdistance2 = rdistancey2+rdistancez2;
194 if (rdistance2>16) return 3;
197 if ((rdistancey2>9.*factor || rdistancez2>9.*factor) && cluster->GetType()==0)
198 return 2; //suspisiouce - will be changed
200 if ((rdistancey2>6.25*factor || rdistancez2>6.25*factor) && cluster->GetType()>0)
201 // strict cut on overlaped cluster
202 return 2; //suspisiouce - will be changed
204 if ( (rdistancey2>1.*factor || rdistancez2>6.25*factor )
205 && cluster->GetType()<0){
215 //_____________________________________________________________________________
216 AliTPCtrackerMI::AliTPCtrackerMI(const AliTPCParam *par):
217 AliTracker(), fkNIS(par->GetNInnerSector()/2), fkNOS(par->GetNOuterSector()/2)
219 //---------------------------------------------------------------------
220 // The main TPC tracker constructor
221 //---------------------------------------------------------------------
222 fInnerSec=new AliTPCSector[fkNIS];
223 fOuterSec=new AliTPCSector[fkNOS];
226 for (i=0; i<fkNIS; i++) fInnerSec[i].Setup(par,0);
227 for (i=0; i<fkNOS; i++) fOuterSec[i].Setup(par,1);
234 Int_t nrowlow = par->GetNRowLow();
235 Int_t nrowup = par->GetNRowUp();
238 for (Int_t i=0;i<nrowlow;i++){
239 fXRow[i] = par->GetPadRowRadiiLow(i);
240 fPadLength[i]= par->GetPadPitchLength(0,i);
241 fYMax[i] = fXRow[i]*TMath::Tan(0.5*par->GetInnerAngle());
245 for (Int_t i=0;i<nrowup;i++){
246 fXRow[i+nrowlow] = par->GetPadRowRadiiUp(i);
247 fPadLength[i+nrowlow] = par->GetPadPitchLength(60,i);
248 fYMax[i+nrowlow] = fXRow[i+nrowlow]*TMath::Tan(0.5*par->GetOuterAngle());
259 fDebugStreamer = new TTreeSRedirector("TPCdebug.root");
261 //________________________________________________________________________
262 AliTPCtrackerMI::AliTPCtrackerMI(const AliTPCtrackerMI &t):
267 //------------------------------------
268 // dummy copy constructor
269 //------------------------------------------------------------------
271 AliTPCtrackerMI & AliTPCtrackerMI::operator=(const AliTPCtrackerMI& /*r*/){
272 //------------------------------
274 //--------------------------------------------------------------
277 //_____________________________________________________________________________
278 AliTPCtrackerMI::~AliTPCtrackerMI() {
279 //------------------------------------------------------------------
280 // TPC tracker destructor
281 //------------------------------------------------------------------
288 if (fDebugStreamer) delete fDebugStreamer;
291 void AliTPCtrackerMI::SetIO()
295 fInput = AliRunLoader::GetTreeR("TPC", kFALSE,AliConfig::GetDefaultEventFolderName());
297 fOutput = AliRunLoader::GetTreeT("TPC", kTRUE,AliConfig::GetDefaultEventFolderName());
299 AliTPCtrack *iotrack= new AliTPCtrack;
300 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
306 void AliTPCtrackerMI::SetIO(TTree * input, TTree * output, AliESD * event)
322 AliTPCtrack *iotrack= new AliTPCtrack;
323 // iotrack->fHelixIn = new TClonesArray("AliHelix");
324 //iotrack->fHelixOut = new TClonesArray("AliHelix");
325 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
328 if (output && (fDebug&2)){
329 //write the full seed information if specified in debug mode
331 fSeedTree = new TTree("Seeds","Seeds");
332 AliTPCseed * vseed = new AliTPCseed;
334 TClonesArray * arrtr = new TClonesArray("AliTPCTrackPoint",160);
335 arrtr->ExpandCreateFast(160);
336 TClonesArray * arre = new TClonesArray("AliTPCExactPoint",160);
338 vseed->fPoints = arrtr;
339 vseed->fEPoints = arre;
340 // vseed->fClusterPoints = arrcl;
341 fSeedTree->Branch("seeds","AliTPCseed",&vseed,32000,99);
344 fTreeDebug = new TTree("trackDebug","trackDebug");
345 TClonesArray * arrd = new TClonesArray("AliTPCTrackPoint2",0);
346 fTreeDebug->Branch("debug",&arrd,32000,99);
354 void AliTPCtrackerMI::FillESD(TObjArray* arr)
358 //fill esds using updated tracks
360 // write tracks to the event
361 // store index of the track
362 Int_t nseed=arr->GetEntriesFast();
363 //FindKinks(arr,fEvent);
364 for (Int_t i=0; i<nseed; i++) {
365 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
368 // pt->PropagateTo(fParam->GetInnerRadiusLow());
369 if (pt->GetKinkIndex(0)<=0){ //don't propagate daughter tracks
370 pt->PropagateTo(fParam->GetInnerRadiusLow());
373 if (( pt->GetPoints()[2]- pt->GetPoints()[0])>5 && pt->GetPoints()[3]>0.8){
375 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
376 iotrack.SetTPCPoints(pt->GetPoints());
377 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
378 iotrack.SetV0Indexes(pt->GetV0Indexes());
379 // iotrack.SetTPCpid(pt->fTPCr);
380 //iotrack.SetTPCindex(i);
381 fEvent->AddTrack(&iotrack);
385 if ( (pt->GetNumberOfClusters()>70)&& (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.55) {
387 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
388 iotrack.SetTPCPoints(pt->GetPoints());
389 //iotrack.SetTPCindex(i);
390 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
391 iotrack.SetV0Indexes(pt->GetV0Indexes());
392 // iotrack.SetTPCpid(pt->fTPCr);
393 fEvent->AddTrack(&iotrack);
397 // short tracks - maybe decays
399 if ( (pt->GetNumberOfClusters()>30) && (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.70) {
400 Int_t found,foundable,shared;
401 pt->GetClusterStatistic(0,60,found, foundable,shared,kFALSE);
402 if ( (found>20) && (pt->fNShared/float(pt->GetNumberOfClusters())<0.2)){
404 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
405 //iotrack.SetTPCindex(i);
406 iotrack.SetTPCPoints(pt->GetPoints());
407 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
408 iotrack.SetV0Indexes(pt->GetV0Indexes());
409 //iotrack.SetTPCpid(pt->fTPCr);
410 fEvent->AddTrack(&iotrack);
415 if ( (pt->GetNumberOfClusters()>20) && (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.8) {
416 Int_t found,foundable,shared;
417 pt->GetClusterStatistic(0,60,found, foundable,shared,kFALSE);
418 if (found<20) continue;
419 if (pt->fNShared/float(pt->GetNumberOfClusters())>0.2) continue;
422 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
423 iotrack.SetTPCPoints(pt->GetPoints());
424 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
425 iotrack.SetV0Indexes(pt->GetV0Indexes());
426 //iotrack.SetTPCpid(pt->fTPCr);
427 //iotrack.SetTPCindex(i);
428 fEvent->AddTrack(&iotrack);
431 // short tracks - secondaties
433 if ( (pt->GetNumberOfClusters()>30) ) {
434 Int_t found,foundable,shared;
435 pt->GetClusterStatistic(128,158,found, foundable,shared,kFALSE);
436 if ( (found>20) && (pt->fNShared/float(pt->GetNumberOfClusters())<0.2) &&float(found)/float(foundable)>0.8){
438 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
439 iotrack.SetTPCPoints(pt->GetPoints());
440 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
441 iotrack.SetV0Indexes(pt->GetV0Indexes());
442 //iotrack.SetTPCpid(pt->fTPCr);
443 //iotrack.SetTPCindex(i);
444 fEvent->AddTrack(&iotrack);
449 if ( (pt->GetNumberOfClusters()>15)) {
450 Int_t found,foundable,shared;
451 pt->GetClusterStatistic(138,158,found, foundable,shared,kFALSE);
452 if (found<15) continue;
453 if (pt->fNShared/float(pt->GetNumberOfClusters())>0.2) continue;
454 if (float(found)/float(foundable)<0.8) continue;
457 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
458 iotrack.SetTPCPoints(pt->GetPoints());
459 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
460 iotrack.SetV0Indexes(pt->GetV0Indexes());
461 // iotrack.SetTPCpid(pt->fTPCr);
462 //iotrack.SetTPCindex(i);
463 fEvent->AddTrack(&iotrack);
468 printf("Number of filled ESDs-\t%d\n",fEvent->GetNumberOfTracks());
471 void AliTPCtrackerMI::WriteTracks(TTree * tree)
474 // write tracks from seed array to selected tree
478 AliTPCtrack *iotrack= new AliTPCtrack;
479 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
484 void AliTPCtrackerMI::WriteTracks()
487 // write tracks to the given output tree -
488 // output specified with SetIO routine
495 AliTPCtrack *iotrack= 0;
496 Int_t nseed=fSeeds->GetEntriesFast();
497 //for (Int_t i=0; i<nseed; i++) {
498 // iotrack= (AliTPCtrack*)fSeeds->UncheckedAt(i);
499 // if (iotrack) break;
501 //TBranch * br = fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
502 TBranch * br = fOutput->GetBranch("tracks");
503 br->SetAddress(&iotrack);
505 for (Int_t i=0; i<nseed; i++) {
506 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i);
508 AliTPCtrack * track = new AliTPCtrack(*pt);
511 // br->SetAddress(&iotrack);
516 //fOutput->GetDirectory()->cd();
522 //write the full seed information if specified in debug mode
524 AliTPCseed * vseed = new AliTPCseed;
526 TClonesArray * arrtr = new TClonesArray("AliTPCTrackPoint",160);
527 arrtr->ExpandCreateFast(160);
528 //TClonesArray * arrcl = new TClonesArray("AliTPCclusterMI",160);
529 //arrcl->ExpandCreateFast(160);
530 TClonesArray * arre = new TClonesArray("AliTPCExactPoint",160);
532 vseed->fPoints = arrtr;
533 vseed->fEPoints = arre;
534 // vseed->fClusterPoints = arrcl;
535 //TBranch * brseed = seedtree->Branch("seeds","AliTPCseed",&vseed,32000,99);
536 TBranch * brseed = fSeedTree->GetBranch("seeds");
538 Int_t nseed=fSeeds->GetEntriesFast();
540 for (Int_t i=0; i<nseed; i++) {
541 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i);
544 // pt->fClusterPoints = arrcl;
548 brseed->SetAddress(&vseed);
552 // pt->fClusterPoints = 0;
555 if (fTreeDebug) fTreeDebug->Write();
563 Double_t AliTPCtrackerMI::ErrY2(AliTPCseed* seed, AliTPCclusterMI * cl){
566 //seed->SetErrorY2(0.1);
568 //calculate look-up table at the beginning
569 static Bool_t ginit = kFALSE;
570 static Float_t gnoise1,gnoise2,gnoise3;
571 static Float_t ggg1[10000];
572 static Float_t ggg2[10000];
573 static Float_t ggg3[10000];
574 static Float_t glandau1[10000];
575 static Float_t glandau2[10000];
576 static Float_t glandau3[10000];
578 static Float_t gcor01[500];
579 static Float_t gcor02[500];
580 static Float_t gcorp[500];
585 for (Int_t i=1;i<500;i++){
586 Float_t rsigma = float(i)/100.;
587 gcor02[i] = TMath::Max(0.78 +TMath::Exp(7.4*(rsigma-1.2)),0.6);
588 gcor01[i] = TMath::Max(0.72 +TMath::Exp(3.36*(rsigma-1.2)),0.6);
589 gcorp[i] = TMath::Max(TMath::Power((rsigma+0.5),1.5),1.2);
593 for (Int_t i=3;i<10000;i++){
597 Float_t amp = float(i);
598 Float_t padlength =0.75;
599 gnoise1 = 0.0004/padlength;
600 Float_t nel = 0.268*amp;
601 Float_t nprim = 0.155*amp;
602 ggg1[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.001*nel/(padlength*padlength))/nel;
603 glandau1[i] = (2.+0.12*nprim)*0.5* (2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
604 if (glandau1[i]>1) glandau1[i]=1;
605 glandau1[i]*=padlength*padlength/12.;
609 gnoise2 = 0.0004/padlength;
612 ggg2[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
613 glandau2[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
614 if (glandau2[i]>1) glandau2[i]=1;
615 glandau2[i]*=padlength*padlength/12.;
620 gnoise3 = 0.0004/padlength;
623 ggg3[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
624 glandau3[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
625 if (glandau3[i]>1) glandau3[i]=1;
626 glandau3[i]*=padlength*padlength/12.;
634 Int_t amp = int(TMath::Abs(cl->GetQ()));
636 seed->SetErrorY2(1.);
640 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
641 Int_t ctype = cl->GetType();
642 Float_t padlength= GetPadPitchLength(seed->fRow);
643 Double_t angle2 = seed->GetSnp()*seed->GetSnp();
644 angle2 = angle2/(1-angle2);
647 Int_t rsigmay = int(100.*cl->GetSigmaY2()/(seed->fCurrentSigmaY2));
650 if (fSectors==fInnerSec){
652 res = ggg1[amp]*z+glandau1[amp]*angle2;
653 if (ctype==0) res *= gcor01[rsigmay];
656 res*= gcorp[rsigmay];
662 res = ggg2[amp]*z+glandau2[amp]*angle2;
663 if (ctype==0) res *= gcor02[rsigmay];
666 res*= gcorp[rsigmay];
671 res = ggg3[amp]*z+glandau3[amp]*angle2;
672 if (ctype==0) res *= gcor02[rsigmay];
675 res*= gcorp[rsigmay];
682 res*=2.4; // overestimate error 2 times
689 seed->SetErrorY2(res);
697 Double_t AliTPCtrackerMI::ErrZ2(AliTPCseed* seed, AliTPCclusterMI * cl){
700 //seed->SetErrorY2(0.1);
702 //calculate look-up table at the beginning
703 static Bool_t ginit = kFALSE;
704 static Float_t gnoise1,gnoise2,gnoise3;
705 static Float_t ggg1[10000];
706 static Float_t ggg2[10000];
707 static Float_t ggg3[10000];
708 static Float_t glandau1[10000];
709 static Float_t glandau2[10000];
710 static Float_t glandau3[10000];
712 static Float_t gcor01[1000];
713 static Float_t gcor02[1000];
714 static Float_t gcorp[1000];
719 for (Int_t i=1;i<1000;i++){
720 Float_t rsigma = float(i)/100.;
721 gcor02[i] = TMath::Max(0.81 +TMath::Exp(6.8*(rsigma-1.2)),0.6);
722 gcor01[i] = TMath::Max(0.72 +TMath::Exp(2.04*(rsigma-1.2)),0.6);
723 gcorp[i] = TMath::Max(TMath::Power((rsigma+0.5),1.5),1.2);
727 for (Int_t i=3;i<10000;i++){
731 Float_t amp = float(i);
732 Float_t padlength =0.75;
733 gnoise1 = 0.0004/padlength;
734 Float_t nel = 0.268*amp;
735 Float_t nprim = 0.155*amp;
736 ggg1[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.001*nel/(padlength*padlength))/nel;
737 glandau1[i] = (2.+0.12*nprim)*0.5* (2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
738 if (glandau1[i]>1) glandau1[i]=1;
739 glandau1[i]*=padlength*padlength/12.;
743 gnoise2 = 0.0004/padlength;
746 ggg2[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
747 glandau2[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
748 if (glandau2[i]>1) glandau2[i]=1;
749 glandau2[i]*=padlength*padlength/12.;
754 gnoise3 = 0.0004/padlength;
757 ggg3[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
758 glandau3[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
759 if (glandau3[i]>1) glandau3[i]=1;
760 glandau3[i]*=padlength*padlength/12.;
768 Int_t amp = int(TMath::Abs(cl->GetQ()));
770 seed->SetErrorY2(1.);
774 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
775 Int_t ctype = cl->GetType();
776 Float_t padlength= GetPadPitchLength(seed->fRow);
778 Double_t angle2 = seed->GetSnp()*seed->GetSnp();
779 // if (angle2<0.6) angle2 = 0.6;
780 angle2 = seed->GetTgl()*seed->GetTgl()*(1+angle2/(1-angle2));
783 Int_t rsigmaz = int(100.*cl->GetSigmaZ2()/(seed->fCurrentSigmaZ2));
786 if (fSectors==fInnerSec){
788 res = ggg1[amp]*z+glandau1[amp]*angle2;
789 if (ctype==0) res *= gcor01[rsigmaz];
792 res*= gcorp[rsigmaz];
798 res = ggg2[amp]*z+glandau2[amp]*angle2;
799 if (ctype==0) res *= gcor02[rsigmaz];
802 res*= gcorp[rsigmaz];
807 res = ggg3[amp]*z+glandau3[amp]*angle2;
808 if (ctype==0) res *= gcor02[rsigmaz];
811 res*= gcorp[rsigmaz];
820 if ((ctype<0) &&<70){
828 seed->SetErrorZ2(res);
835 Double_t AliTPCtrackerMI::ErrZ2(AliTPCseed* seed, AliTPCclusterMI * cl){
838 //seed->SetErrorZ2(0.1);
842 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
844 Float_t rsigmaz = cl->GetSigmaZ2()/(seed->fCurrentSigmaZ2);
845 Int_t ctype = cl->GetType();
846 Float_t amp = TMath::Abs(cl->GetQ());
851 Float_t landau=2 ; //landau fluctuation part
852 Float_t gg=2; // gg fluctuation part
853 Float_t padlength= GetPadPitchLength(seed->GetX());
855 if (fSectors==fInnerSec){
856 snoise2 = 0.0004/padlength;
859 gg = (2+0.001*nel/(padlength*padlength))/nel;
860 landau = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
861 if (landau>1) landau=1;
864 snoise2 = 0.0004/padlength;
867 gg = (2+0.0008*nel/(padlength*padlength))/nel;
868 landau = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
869 if (landau>1) landau=1;
871 Float_t sdiff = gg*fParam->GetDiffT()*fParam->GetDiffT()*z;
874 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
875 angle2 = TMath::Sqrt((1-angle2));
876 if (angle2<0.6) angle2 = 0.6;
879 Float_t angle = seed->GetTgl()/angle2;
880 Float_t angular = landau*angle*angle*padlength*padlength/12.;
881 Float_t res = sdiff + angular;
884 if ((ctype==0) && (fSectors ==fOuterSec))
885 res *= 0.81 +TMath::Exp(6.8*(rsigmaz-1.2));
887 if ((ctype==0) && (fSectors ==fInnerSec))
888 res *= 0.72 +TMath::Exp(2.04*(rsigmaz-1.2));
892 res*= TMath::Power(rsigmaz+0.5,1.5); //0.31+0.147*ctype;
898 if ((ctype<0) &&<70){
906 seed->SetErrorZ2(res);
914 void AliTPCtrackerMI::RotateToLocal(AliTPCseed *seed)
916 //rotate to track "local coordinata
917 Float_t x = seed->GetX();
918 Float_t y = seed->GetY();
919 Float_t ymax = x*TMath::Tan(0.5*fSectors->GetAlpha());
922 seed->fRelativeSector= (seed->fRelativeSector+1) % fN;
923 if (!seed->Rotate(fSectors->GetAlpha()))
925 } else if (y <-ymax) {
926 seed->fRelativeSector= (seed->fRelativeSector-1+fN) % fN;
927 if (!seed->Rotate(-fSectors->GetAlpha()))
935 //_____________________________________________________________________________
936 Double_t AliTPCtrackerMI::F1old(Double_t x1,Double_t y1,
937 Double_t x2,Double_t y2,
938 Double_t x3,Double_t y3)
940 //-----------------------------------------------------------------
941 // Initial approximation of the track curvature
942 //-----------------------------------------------------------------
943 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
944 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
945 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
946 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
947 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
949 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
950 if ( xr*xr+yr*yr<=0.00000000000001) return 100;
951 return -xr*yr/sqrt(xr*xr+yr*yr);
956 //_____________________________________________________________________________
957 Double_t AliTPCtrackerMI::F1(Double_t x1,Double_t y1,
958 Double_t x2,Double_t y2,
959 Double_t x3,Double_t y3)
961 //-----------------------------------------------------------------
962 // Initial approximation of the track curvature
963 //-----------------------------------------------------------------
969 Double_t det = x3*y2-x2*y3;
974 Double_t u = 0.5* (x2*(x2-x3)+y2*(y2-y3))/det;
975 Double_t x0 = x3*0.5-y3*u;
976 Double_t y0 = y3*0.5+x3*u;
977 Double_t c2 = 1/TMath::Sqrt(x0*x0+y0*y0);
983 Double_t AliTPCtrackerMI::F2(Double_t x1,Double_t y1,
984 Double_t x2,Double_t y2,
985 Double_t x3,Double_t y3)
987 //-----------------------------------------------------------------
988 // Initial approximation of the track curvature
989 //-----------------------------------------------------------------
995 Double_t det = x3*y2-x2*y3;
1000 Double_t u = 0.5* (x2*(x2-x3)+y2*(y2-y3))/det;
1001 Double_t x0 = x3*0.5-y3*u;
1002 Double_t y0 = y3*0.5+x3*u;
1003 Double_t c2 = 1/TMath::Sqrt(x0*x0+y0*y0);
1012 //_____________________________________________________________________________
1013 Double_t AliTPCtrackerMI::F2old(Double_t x1,Double_t y1,
1014 Double_t x2,Double_t y2,
1015 Double_t x3,Double_t y3)
1017 //-----------------------------------------------------------------
1018 // Initial approximation of the track curvature times center of curvature
1019 //-----------------------------------------------------------------
1020 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
1021 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
1022 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
1023 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
1024 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
1026 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
1028 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
1031 //_____________________________________________________________________________
1032 Double_t AliTPCtrackerMI::F3(Double_t x1,Double_t y1,
1033 Double_t x2,Double_t y2,
1034 Double_t z1,Double_t z2)
1036 //-----------------------------------------------------------------
1037 // Initial approximation of the tangent of the track dip angle
1038 //-----------------------------------------------------------------
1039 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
1043 Double_t AliTPCtrackerMI::F3n(Double_t x1,Double_t y1,
1044 Double_t x2,Double_t y2,
1045 Double_t z1,Double_t z2, Double_t c)
1047 //-----------------------------------------------------------------
1048 // Initial approximation of the tangent of the track dip angle
1049 //-----------------------------------------------------------------
1053 //angle1 = (z1-z2)*c/(TMath::ASin(c*x1-ni)-TMath::ASin(c*x2-ni));
1055 Double_t d = TMath::Sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
1056 if (TMath::Abs(d*c*0.5)>1) return 0;
1057 // Double_t angle2 = TMath::ASin(d*c*0.5);
1058 // Double_t angle2 = AliTPCFastMath::FastAsin(d*c*0.5);
1059 Double_t angle2 = (d*c*0.5>0.1)? TMath::ASin(d*c*0.5): AliTPCFastMath::FastAsin(d*c*0.5);
1061 angle2 = (z1-z2)*c/(angle2*2.);
1065 Bool_t AliTPCtrackerMI::GetProlongation(Double_t x1, Double_t x2, Double_t x[5], Double_t &y, Double_t &z)
1066 {//-----------------------------------------------------------------
1067 // This function find proloncation of a track to a reference plane x=x2.
1068 //-----------------------------------------------------------------
1072 if (TMath::Abs(x[4]*x1 - x[2]) >= 0.999) {
1076 Double_t c1=x[4]*x1 - x[2], r1=sqrt(1.- c1*c1);
1077 Double_t c2=x[4]*x2 - x[2], r2=sqrt(1.- c2*c2);
1081 Double_t dy = dx*(c1+c2)/(r1+r2);
1084 Double_t delta = x[4]*dx*(c1+c2)/(c1*r2 + c2*r1);
1086 if (TMath::Abs(delta)>0.01){
1087 dz = x[3]*TMath::ASin(delta)/x[4];
1089 dz = x[3]*AliTPCFastMath::FastAsin(delta)/x[4];
1092 //dz = x[3]*AliTPCFastMath::FastAsin(delta)/x[4];
1100 Int_t AliTPCtrackerMI::LoadClusters (TTree *tree)
1105 return LoadClusters();
1108 Int_t AliTPCtrackerMI::LoadClusters()
1111 // load clusters to the memory
1112 AliTPCClustersRow *clrow= new AliTPCClustersRow;
1113 clrow->SetClass("AliTPCclusterMI");
1115 clrow->GetArray()->ExpandCreateFast(10000);
1117 // TTree * tree = fClustersArray.GetTree();
1119 TTree * tree = fInput;
1120 TBranch * br = tree->GetBranch("Segment");
1121 br->SetAddress(&clrow);
1123 Int_t j=Int_t(tree->GetEntries());
1124 for (Int_t i=0; i<j; i++) {
1128 fParam->AdjustSectorRow(clrow->GetID(),sec,row);
1129 for (Int_t icl=0; icl<clrow->GetArray()->GetEntriesFast(); icl++){
1130 Transform((AliCluster*)(clrow->GetArray()->At(icl)));
1133 AliTPCRow * tpcrow=0;
1136 tpcrow = &(fInnerSec[sec%fkNIS][row]);
1140 tpcrow = &(fOuterSec[(sec-fkNIS*2)%fkNOS][row]);
1141 left = (sec-fkNIS*2)/fkNOS;
1144 tpcrow->fN1 = clrow->GetArray()->GetEntriesFast();
1145 tpcrow->fClusters1 = new AliTPCclusterMI[tpcrow->fN1];
1146 for (Int_t i=0;i<tpcrow->fN1;i++)
1147 tpcrow->fClusters1[i] = *(AliTPCclusterMI*)(clrow->GetArray()->At(i));
1150 tpcrow->fN2 = clrow->GetArray()->GetEntriesFast();
1151 tpcrow->fClusters2 = new AliTPCclusterMI[tpcrow->fN2];
1152 for (Int_t i=0;i<tpcrow->fN2;i++)
1153 tpcrow->fClusters2[i] = *(AliTPCclusterMI*)(clrow->GetArray()->At(i));
1164 void AliTPCtrackerMI::UnloadClusters()
1167 // unload clusters from the memory
1169 Int_t nrows = fOuterSec->GetNRows();
1170 for (Int_t sec = 0;sec<fkNOS;sec++)
1171 for (Int_t row = 0;row<nrows;row++){
1172 AliTPCRow* tpcrow = &(fOuterSec[sec%fkNOS][row]);
1174 // if (tpcrow->fClusters1) delete []tpcrow->fClusters1;
1175 // if (tpcrow->fClusters2) delete []tpcrow->fClusters2;
1177 tpcrow->ResetClusters();
1180 nrows = fInnerSec->GetNRows();
1181 for (Int_t sec = 0;sec<fkNIS;sec++)
1182 for (Int_t row = 0;row<nrows;row++){
1183 AliTPCRow* tpcrow = &(fInnerSec[sec%fkNIS][row]);
1185 // if (tpcrow->fClusters1) delete []tpcrow->fClusters1;
1186 //if (tpcrow->fClusters2) delete []tpcrow->fClusters2;
1188 tpcrow->ResetClusters();
1194 void AliTPCtrackerMI::Transform(AliCluster * cluster){
1198 //if (!fParam->IsGeoRead()) fParam->ReadGeoMatrices();
1199 TGeoHMatrix *mat = fParam->GetClusterMatrix(cluster->GetDetector());
1200 Double_t pos[3]= {cluster->GetX(),cluster->GetY(),cluster->GetZ()};
1202 mat->LocalToMaster(pos,posC);
1203 cluster->SetX(posC[0]);
1204 cluster->SetY(posC[1]);
1205 cluster->SetZ(posC[2]);
1208 //_____________________________________________________________________________
1209 Int_t AliTPCtrackerMI::LoadOuterSectors() {
1210 //-----------------------------------------------------------------
1211 // This function fills outer TPC sectors with clusters.
1212 //-----------------------------------------------------------------
1213 Int_t nrows = fOuterSec->GetNRows();
1215 for (Int_t sec = 0;sec<fkNOS;sec++)
1216 for (Int_t row = 0;row<nrows;row++){
1217 AliTPCRow* tpcrow = &(fOuterSec[sec%fkNOS][row]);
1218 Int_t sec2 = sec+2*fkNIS;
1220 Int_t ncl = tpcrow->fN1;
1222 AliTPCclusterMI *c= &(tpcrow->fClusters1[ncl]);
1223 index=(((sec2<<8)+row)<<16)+ncl;
1224 tpcrow->InsertCluster(c,index);
1229 AliTPCclusterMI *c= &(tpcrow->fClusters2[ncl]);
1230 index=((((sec2+fkNOS)<<8)+row)<<16)+ncl;
1231 tpcrow->InsertCluster(c,index);
1234 // write indexes for fast acces
1236 for (Int_t i=0;i<510;i++)
1237 tpcrow->fFastCluster[i]=-1;
1238 for (Int_t i=0;i<tpcrow->GetN();i++){
1239 Int_t zi = Int_t((*tpcrow)[i]->GetZ()+255.);
1240 tpcrow->fFastCluster[zi]=i; // write index
1243 for (Int_t i=0;i<510;i++){
1244 if (tpcrow->fFastCluster[i]<0)
1245 tpcrow->fFastCluster[i] = last;
1247 last = tpcrow->fFastCluster[i];
1256 //_____________________________________________________________________________
1257 Int_t AliTPCtrackerMI::LoadInnerSectors() {
1258 //-----------------------------------------------------------------
1259 // This function fills inner TPC sectors with clusters.
1260 //-----------------------------------------------------------------
1261 Int_t nrows = fInnerSec->GetNRows();
1263 for (Int_t sec = 0;sec<fkNIS;sec++)
1264 for (Int_t row = 0;row<nrows;row++){
1265 AliTPCRow* tpcrow = &(fInnerSec[sec%fkNIS][row]);
1268 Int_t ncl = tpcrow->fN1;
1270 AliTPCclusterMI *c= &(tpcrow->fClusters1[ncl]);
1271 index=(((sec<<8)+row)<<16)+ncl;
1272 tpcrow->InsertCluster(c,index);
1277 AliTPCclusterMI *c= &(tpcrow->fClusters2[ncl]);
1278 index=((((sec+fkNIS)<<8)+row)<<16)+ncl;
1279 tpcrow->InsertCluster(c,index);
1282 // write indexes for fast acces
1284 for (Int_t i=0;i<510;i++)
1285 tpcrow->fFastCluster[i]=-1;
1286 for (Int_t i=0;i<tpcrow->GetN();i++){
1287 Int_t zi = Int_t((*tpcrow)[i]->GetZ()+255.);
1288 tpcrow->fFastCluster[zi]=i; // write index
1291 for (Int_t i=0;i<510;i++){
1292 if (tpcrow->fFastCluster[i]<0)
1293 tpcrow->fFastCluster[i] = last;
1295 last = tpcrow->fFastCluster[i];
1307 //_________________________________________________________________________
1308 AliTPCclusterMI *AliTPCtrackerMI::GetClusterMI(Int_t index) const {
1309 //--------------------------------------------------------------------
1310 // Return pointer to a given cluster
1311 //--------------------------------------------------------------------
1312 Int_t sec=(index&0xff000000)>>24;
1313 Int_t row=(index&0x00ff0000)>>16;
1314 Int_t ncl=(index&0x00007fff)>>00;
1316 const AliTPCRow * tpcrow=0;
1317 AliTPCclusterMI * clrow =0;
1319 if (sec<0 || sec>=fkNIS*4) {
1320 AliWarning(Form("Wrong sector %d",sec));
1325 tpcrow = &(fInnerSec[sec%fkNIS][row]);
1326 if (tpcrow==0) return 0;
1329 if (tpcrow->fN1<=ncl) return 0;
1330 clrow = tpcrow->fClusters1;
1333 if (tpcrow->fN2<=ncl) return 0;
1334 clrow = tpcrow->fClusters2;
1338 tpcrow = &(fOuterSec[(sec-fkNIS*2)%fkNOS][row]);
1339 if (tpcrow==0) return 0;
1341 if (sec-2*fkNIS<fkNOS) {
1342 if (tpcrow->fN1<=ncl) return 0;
1343 clrow = tpcrow->fClusters1;
1346 if (tpcrow->fN2<=ncl) return 0;
1347 clrow = tpcrow->fClusters2;
1351 return &(clrow[ncl]);
1357 Int_t AliTPCtrackerMI::FollowToNext(AliTPCseed& t, Int_t nr) {
1358 //-----------------------------------------------------------------
1359 // This function tries to find a track prolongation to next pad row
1360 //-----------------------------------------------------------------
1362 Double_t x= GetXrow(nr), ymax=GetMaxY(nr);
1363 AliTPCclusterMI *cl=0;
1364 Int_t tpcindex= t.GetClusterIndex2(nr);
1366 // update current shape info every 5 pad-row
1367 // if ( (nr%5==0) || t.GetNumberOfClusters()<2 || (t.fCurrentSigmaY2<0.0001) ){
1371 if (fIteration>0 && tpcindex>=-1){ //if we have already clusters
1373 if (tpcindex==-1) return 0; //track in dead zone
1375 cl = t.fClusterPointer[nr];
1376 if ( (cl==0) ) cl = GetClusterMI(tpcindex);
1377 t.fCurrentClusterIndex1 = tpcindex;
1380 Int_t relativesector = ((tpcindex&0xff000000)>>24)%18; // if previously accepted cluster in different sector
1381 Float_t angle = relativesector*fSectors->GetAlpha()+fSectors->GetAlphaShift();
1383 if (angle<-TMath::Pi()) angle += 2*TMath::Pi();
1384 if (angle>=TMath::Pi()) angle -= 2*TMath::Pi();
1386 if (TMath::Abs(angle-t.GetAlpha())>0.001){
1387 Double_t rotation = angle-t.GetAlpha();
1388 t.fRelativeSector= relativesector;
1389 if (!t.Rotate(rotation)) return 0;
1391 if (!t.PropagateTo(x)) return 0;
1393 t.fCurrentCluster = cl;
1395 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1396 if ((tpcindex&0x8000)==0) accept =0;
1398 //if founded cluster is acceptible
1399 if (cl->IsUsed(11)) { // id cluster is shared inrease uncertainty
1406 UpdateTrack(&t,accept);
1411 if (TMath::Abs(t.GetSnp())>AliTPCReconstructor::GetMaxSnpTracker()) return 0; // cut on angle
1413 // not look for new cluster during refitting
1419 // if (TMath::Abs(t.GetSnp())>0.95 || TMath::Abs(x*t.GetC()-t.GetEta())>0.95) return 0;// patch 28 fev 06
1420 Double_t y=t.GetYat(x);
1421 if (TMath::Abs(y)>ymax){
1423 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1424 if (!t.Rotate(fSectors->GetAlpha()))
1426 } else if (y <-ymax) {
1427 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1428 if (!t.Rotate(-fSectors->GetAlpha()))
1434 if (!t.PropagateTo(x)) {
1435 if (fIteration==0) t.fRemoval = 10;
1439 Double_t z=t.GetZ();
1441 const AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1442 if ( (t.GetSigmaY2()<0) || t.GetSigmaZ2()<0) return 0;
1444 Double_t roadz = 1.;
1446 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1448 t.SetClusterIndex2(nr,-1);
1453 if (TMath::Abs(z)<(AliTPCReconstructor::GetCtgRange()*x+10) && TMath::Abs(z)<fParam->GetZLength() && (TMath::Abs(t.GetSnp())<AliTPCReconstructor::GetMaxSnpTracker()))
1460 // cl = krow.FindNearest2(y+10.,z,roady,roadz,index);
1461 cl = krow.FindNearest2(y,z,roady,roadz,index);
1462 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1465 t.fCurrentCluster = cl;
1467 if (fIteration==2&&cl->IsUsed(10)) return 0;
1468 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1469 if (fIteration==2&&cl->IsUsed(11)) {
1476 if (t.fCurrentCluster->IsUsed(10)){
1481 if (t.fNShared>0.7*t.GetNumberOfClusters()) {
1487 if (accept<3) UpdateTrack(&t,accept);
1490 if ( fIteration==0 && t.fNFoundable*0.5 > t.GetNumberOfClusters()) t.fRemoval=10;
1496 Int_t AliTPCtrackerMI::FollowToNextFast(AliTPCseed& t, Int_t nr) {
1497 //-----------------------------------------------------------------
1498 // This function tries to find a track prolongation to next pad row
1499 //-----------------------------------------------------------------
1501 Double_t x= GetXrow(nr), ymax=GetMaxY(nr);
1503 if (!t.GetProlongation(x,y,z)) {
1509 if (TMath::Abs(y)>ymax){
1512 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1513 if (!t.Rotate(fSectors->GetAlpha()))
1515 } else if (y <-ymax) {
1516 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1517 if (!t.Rotate(-fSectors->GetAlpha()))
1520 if (!t.PropagateTo(x)) {
1523 t.GetProlongation(x,y,z);
1526 // update current shape info every 3 pad-row
1527 if ( (nr%6==0) || t.GetNumberOfClusters()<2 || (t.fCurrentSigmaY2<0.0001) ){
1528 // t.fCurrentSigmaY = GetSigmaY(&t);
1529 //t.fCurrentSigmaZ = GetSigmaZ(&t);
1533 AliTPCclusterMI *cl=0;
1538 const AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1539 if ( (t.GetSigmaY2()<0) || t.GetSigmaZ2()<0) return 0;
1541 Double_t roadz = 1.;
1544 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1546 t.SetClusterIndex2(row,-1);
1551 if (TMath::Abs(z)>(AliTPCReconstructor::GetCtgRange()*x+10)) t.SetClusterIndex2(row,-1);
1555 if ((cl==0)&&(krow)) {
1556 // cl = krow.FindNearest2(y+10,z,roady,roadz,index);
1557 cl = krow.FindNearest2(y,z,roady,roadz,index);
1559 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1563 t.fCurrentCluster = cl;
1564 // Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1566 t.SetClusterIndex2(row,index);
1567 t.fClusterPointer[row] = cl;
1575 //_________________________________________________________________________
1576 Bool_t AliTPCtrackerMI::GetTrackPoint(Int_t index, AliTrackPoint &p ) const
1578 // Get track space point by index
1579 // return false in case the cluster doesn't exist
1580 AliTPCclusterMI *cl = GetClusterMI(index);
1581 if (!cl) return kFALSE;
1582 Int_t sector = (index&0xff000000)>>24;
1583 // Int_t row = (index&0x00ff0000)>>16;
1585 // xyz[0] = fParam->GetPadRowRadii(sector,row);
1586 xyz[0] = cl->GetX();
1587 xyz[1] = cl->GetY();
1588 xyz[2] = cl->GetZ();
1590 fParam->AdjustCosSin(sector,cos,sin);
1591 Float_t x = cos*xyz[0]-sin*xyz[1];
1592 Float_t y = cos*xyz[1]+sin*xyz[0];
1594 Float_t sigmaY2 = 0.027*cl->GetSigmaY2();
1595 if (sector < fParam->GetNInnerSector()) sigmaY2 *= 2.07;
1596 Float_t sigmaZ2 = 0.066*cl->GetSigmaZ2();
1597 if (sector < fParam->GetNInnerSector()) sigmaZ2 *= 1.77;
1598 cov[0] = sin*sin*sigmaY2;
1599 cov[1] = -sin*cos*sigmaY2;
1601 cov[3] = cos*cos*sigmaY2;
1604 p.SetXYZ(x,y,xyz[2],cov);
1605 AliAlignObj::ELayerID iLayer;
1607 if (sector < fParam->GetNInnerSector()) {
1608 iLayer = AliAlignObj::kTPC1;
1612 iLayer = AliAlignObj::kTPC2;
1613 idet = sector - fParam->GetNInnerSector();
1615 UShort_t volid = AliAlignObj::LayerToVolUID(iLayer,idet);
1616 p.SetVolumeID(volid);
1622 Int_t AliTPCtrackerMI::UpdateClusters(AliTPCseed& t, Int_t nr) {
1623 //-----------------------------------------------------------------
1624 // This function tries to find a track prolongation to next pad row
1625 //-----------------------------------------------------------------
1626 t.fCurrentCluster = 0;
1627 t.fCurrentClusterIndex1 = 0;
1629 Double_t xt=t.GetX();
1630 Int_t row = GetRowNumber(xt)-1;
1631 Double_t ymax= GetMaxY(nr);
1633 if (row < nr) return 1; // don't prolongate if not information until now -
1634 // if (TMath::Abs(t.GetSnp())>0.9 && t.GetNumberOfClusters()>40. && fIteration!=2) {
1636 // return 0; // not prolongate strongly inclined tracks
1638 // if (TMath::Abs(t.GetSnp())>0.95) {
1640 // return 0; // not prolongate strongly inclined tracks
1641 // }// patch 28 fev 06
1643 Double_t x= GetXrow(nr);
1645 //t.PropagateTo(x+0.02);
1646 //t.PropagateTo(x+0.01);
1647 if (!t.PropagateTo(x)){
1654 if (TMath::Abs(y)>ymax){
1656 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1657 if (!t.Rotate(fSectors->GetAlpha()))
1659 } else if (y <-ymax) {
1660 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1661 if (!t.Rotate(-fSectors->GetAlpha()))
1664 // if (!t.PropagateTo(x)){
1671 if (TMath::Abs(t.GetSnp())>AliTPCReconstructor::GetMaxSnpTracker()) return 0;
1672 AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1674 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1676 t.SetClusterIndex2(nr,-1);
1681 if (TMath::Abs(t.GetZ())<(AliTPCReconstructor::GetCtgRange()*t.GetX()+10) && (TMath::Abs(t.GetSnp())<AliTPCReconstructor::GetMaxSnpTracker())) t.fNFoundable++;
1687 if ( (nr%6==0) || t.GetNumberOfClusters()<2){
1688 // t.fCurrentSigmaY = GetSigmaY(&t);
1689 //t.fCurrentSigmaZ = GetSigmaZ(&t);
1693 AliTPCclusterMI *cl=0;
1696 Double_t roady = 1.;
1697 Double_t roadz = 1.;
1701 index = t.GetClusterIndex2(nr);
1702 if ( (index>0) && (index&0x8000)==0){
1703 cl = t.fClusterPointer[nr];
1704 if ( (cl==0) && (index>0)) cl = GetClusterMI(index);
1705 t.fCurrentClusterIndex1 = index;
1707 t.fCurrentCluster = cl;
1713 // if (index<0) return 0;
1714 UInt_t uindex = TMath::Abs(index);
1717 //cl = krow.FindNearest2(y+10,z,roady,roadz,uindex);
1718 cl = krow.FindNearest2(y,z,roady,roadz,uindex);
1721 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(uindex);
1722 t.fCurrentCluster = cl;
1728 Int_t AliTPCtrackerMI::FollowToNextCluster(AliTPCseed & t, Int_t nr) {
1729 //-----------------------------------------------------------------
1730 // This function tries to find a track prolongation to next pad row
1731 //-----------------------------------------------------------------
1733 //update error according neighborhoud
1735 if (t.fCurrentCluster) {
1737 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1739 if (t.fCurrentCluster->IsUsed(10)){
1745 if (t.fNShared>0.7*t.GetNumberOfClusters()) {
1750 if (fIteration>0) accept = 0;
1751 if (accept<3) UpdateTrack(&t,accept);
1755 if ( ( (t.GetSigmaY2()+t.GetSigmaZ2())>0.16)&& t.GetNumberOfClusters()>18) t.fRemoval=10;
1756 if ( t.GetChi2()/t.GetNumberOfClusters()>6 &&t.GetNumberOfClusters()>18) t.fRemoval=10;
1758 if (( (t.fNFoundable*0.5 > t.GetNumberOfClusters()) || t.fNoCluster>15)) t.fRemoval=10;
1766 //_____________________________________________________________________________
1767 Int_t AliTPCtrackerMI::FollowProlongation(AliTPCseed& t, Int_t rf, Int_t step) {
1768 //-----------------------------------------------------------------
1769 // This function tries to find a track prolongation.
1770 //-----------------------------------------------------------------
1771 Double_t xt=t.GetX();
1773 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1774 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1775 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1777 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1779 Int_t first = GetRowNumber(xt)-1;
1780 for (Int_t nr= first; nr>=rf; nr-=step) {
1782 if (t.GetKinkIndexes()[0]>0){
1783 for (Int_t i=0;i<3;i++){
1784 Int_t index = t.GetKinkIndexes()[i];
1785 if (index==0) break;
1786 if (index<0) continue;
1788 AliESDkink * kink = fEvent->GetKink(index-1);
1790 printf("PROBLEM\n");
1793 Int_t kinkrow = kink->GetTPCRow0()+2+Int_t(0.5/(0.05+kink->GetAngle(2)));
1795 AliExternalTrackParam paramd(t);
1796 kink->SetDaughter(paramd);
1797 kink->SetStatus(2,5);
1804 if (nr==80) t.UpdateReference();
1805 if (nr<fInnerSec->GetNRows())
1806 fSectors = fInnerSec;
1808 fSectors = fOuterSec;
1809 if (FollowToNext(t,nr)==0)
1818 //_____________________________________________________________________________
1819 Int_t AliTPCtrackerMI::FollowProlongationFast(AliTPCseed& t, Int_t rf, Int_t step) {
1820 //-----------------------------------------------------------------
1821 // This function tries to find a track prolongation.
1822 //-----------------------------------------------------------------
1823 Double_t xt=t.GetX();
1825 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1826 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1827 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1828 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1830 for (Int_t nr=GetRowNumber(xt)-1; nr>=rf; nr-=step) {
1832 if (FollowToNextFast(t,nr)==0)
1833 if (!t.IsActive()) return 0;
1843 Int_t AliTPCtrackerMI::FollowBackProlongation(AliTPCseed& t, Int_t rf) {
1844 //-----------------------------------------------------------------
1845 // This function tries to find a track prolongation.
1846 //-----------------------------------------------------------------
1848 Double_t xt=t.GetX();
1849 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1850 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1851 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1852 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1854 Int_t first = t.fFirstPoint;
1855 if (first<GetRowNumber(xt)+1) first = GetRowNumber(xt)+1;
1857 if (first<0) first=0;
1858 for (Int_t nr=first; nr<=rf; nr++) {
1859 // if ( (TMath::Abs(t.GetSnp())>0.95)) break;//patch 28 fev 06
1860 if (t.GetKinkIndexes()[0]<0){
1861 for (Int_t i=0;i<3;i++){
1862 Int_t index = t.GetKinkIndexes()[i];
1863 if (index==0) break;
1864 if (index>0) continue;
1865 index = TMath::Abs(index);
1866 AliESDkink * kink = fEvent->GetKink(index-1);
1868 printf("PROBLEM\n");
1871 Int_t kinkrow = kink->GetTPCRow0()-2-Int_t(0.5/(0.05+kink->GetAngle(2)));
1873 AliExternalTrackParam paramm(t);
1874 kink->SetMother(paramm);
1875 kink->SetStatus(2,1);
1882 if (nr<fInnerSec->GetNRows())
1883 fSectors = fInnerSec;
1885 fSectors = fOuterSec;
1895 Float_t AliTPCtrackerMI::OverlapFactor(AliTPCseed * s1, AliTPCseed * s2, Int_t &sum1, Int_t & sum2)
1903 Float_t dz2 =(s1->GetZ() - s2->GetZ());
1906 Float_t dy2 =TMath::Abs((s1->GetY() - s2->GetY()));
1908 Float_t distance = TMath::Sqrt(dz2+dy2);
1909 if (distance>4.) return 0; // if there are far away - not overlap - to reduce combinatorics
1912 Int_t firstpoint = TMath::Min(s1->fFirstPoint,s2->fFirstPoint);
1913 Int_t lastpoint = TMath::Max(s1->fLastPoint,s2->fLastPoint);
1918 if (firstpoint>lastpoint) {
1919 firstpoint =lastpoint;
1924 for (Int_t i=firstpoint-1;i<lastpoint+1;i++){
1925 if (s1->GetClusterIndex2(i)>0) sum1++;
1926 if (s2->GetClusterIndex2(i)>0) sum2++;
1927 if (s1->GetClusterIndex2(i)==s2->GetClusterIndex2(i) && s1->GetClusterIndex2(i)>0) {
1931 if (sum<5) return 0;
1933 Float_t summin = TMath::Min(sum1+1,sum2+1);
1934 Float_t ratio = (sum+1)/Float_t(summin);
1938 void AliTPCtrackerMI::SignShared(AliTPCseed * s1, AliTPCseed * s2)
1942 if (TMath::Abs(s1->GetC()-s2->GetC())>0.004) return;
1943 if (TMath::Abs(s1->GetTgl()-s2->GetTgl())>0.6) return;
1945 Float_t dz2 =(s1->GetZ() - s2->GetZ());
1947 Float_t dy2 =(s1->GetY() - s2->GetY());
1949 Float_t distance = dz2+dy2;
1950 if (distance>325.) return ; // if there are far away - not overlap - to reduce combinatorics
1955 Int_t firstpoint = TMath::Max(s1->fFirstPoint,s2->fFirstPoint);
1956 Int_t lastpoint = TMath::Min(s1->fLastPoint,s2->fLastPoint);
1958 if (firstpoint>=lastpoint-5) return;;
1960 for (Int_t i=firstpoint;i<lastpoint;i++){
1961 // if ( (s1->GetClusterIndex2(i)&0xFFFF8FFF)==(s2->GetClusterIndex2(i)&0xFFFF8FFF) && s1->GetClusterIndex2(i)>0) {
1962 if ( (s1->GetClusterIndex2(i))==(s2->GetClusterIndex2(i)) && s1->GetClusterIndex2(i)>0) {
1969 for (Int_t i=firstpoint;i<lastpoint;i++){
1970 // if ( (s1->GetClusterIndex2(i)&0xFFFF8FFF)==(s2->GetClusterIndex2(i)&0xFFFF8FFF) && s1->GetClusterIndex2(i)>0) {
1971 if ( (s1->GetClusterIndex2(i))==(s2->GetClusterIndex2(i)) && s1->GetClusterIndex2(i)>0) {
1972 AliTPCTrackerPoint *p1 = s1->GetTrackPoint(i);
1973 AliTPCTrackerPoint *p2 = s2->GetTrackPoint(i);;
1974 if (s1->IsActive()&&s2->IsActive()){
1975 p1->fIsShared = kTRUE;
1976 p2->fIsShared = kTRUE;
1983 for (Int_t i=0;i<4;i++){
1984 if (s1->fOverlapLabels[3*i]==0){
1985 s1->fOverlapLabels[3*i] = s2->GetLabel();
1986 s1->fOverlapLabels[3*i+1] = sumshared;
1987 s1->fOverlapLabels[3*i+2] = s2->GetUniqueID();
1991 for (Int_t i=0;i<4;i++){
1992 if (s2->fOverlapLabels[3*i]==0){
1993 s2->fOverlapLabels[3*i] = s1->GetLabel();
1994 s2->fOverlapLabels[3*i+1] = sumshared;
1995 s2->fOverlapLabels[3*i+2] = s1->GetUniqueID();
2003 void AliTPCtrackerMI::SignShared(TObjArray * arr)
2006 //sort trackss according sectors
2008 for (Int_t i=0; i<arr->GetEntriesFast(); i++) {
2009 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2011 //if (pt) RotateToLocal(pt);
2015 arr->Sort(); // sorting according z
2016 arr->Expand(arr->GetEntries());
2019 Int_t nseed=arr->GetEntriesFast();
2020 for (Int_t i=0; i<nseed; i++) {
2021 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2023 for (Int_t j=0;j<=12;j++){
2024 pt->fOverlapLabels[j] =0;
2027 for (Int_t i=0; i<nseed; i++) {
2028 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2030 if (pt->fRemoval>10) continue;
2031 for (Int_t j=i+1; j<nseed; j++){
2032 AliTPCseed *pt2=(AliTPCseed*)arr->UncheckedAt(j);
2034 if (pt2->fRemoval<=10) {
2035 if ( TMath::Abs(pt->fRelativeSector-pt2->fRelativeSector)>0) break;
2042 void AliTPCtrackerMI::RemoveDouble(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t removalindex)
2045 //sort trackss according sectors
2048 Info("RemoveDouble","Number of tracks before double removal- %d\n",arr->GetEntries());
2051 for (Int_t i=0; i<arr->GetEntriesFast(); i++) {
2052 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2057 arr->Sort(); // sorting according z
2058 arr->Expand(arr->GetEntries());
2060 //reset overlap labels
2062 Int_t nseed=arr->GetEntriesFast();
2063 for (Int_t i=0; i<nseed; i++) {
2064 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2067 for (Int_t j=0;j<=12;j++){
2068 pt->fOverlapLabels[j] =0;
2072 //sign shared tracks
2073 for (Int_t i=0; i<nseed; i++) {
2074 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2076 if (pt->fRemoval>10) continue;
2077 Float_t deltac = pt->GetC()*0.1;
2078 for (Int_t j=i+1; j<nseed; j++){
2079 AliTPCseed *pt2=(AliTPCseed*)arr->UncheckedAt(j);
2081 if (pt2->fRemoval<=10) {
2082 if ( TMath::Abs(pt->fRelativeSector-pt2->fRelativeSector)>0) break;
2083 if (TMath::Abs(pt->GetC() -pt2->GetC())>deltac) continue;
2084 if (TMath::Abs(pt->GetTgl()-pt2->GetTgl())>0.05) continue;
2091 // remove highly shared tracks
2092 for (Int_t i=0; i<nseed; i++) {
2093 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2095 if (pt->fRemoval>10) continue;
2098 for (Int_t j=0;j<4;j++){
2099 sumshared = pt->fOverlapLabels[j*3+1];
2101 Float_t factor = factor1;
2102 if (pt->fRemoval>0) factor = factor2;
2103 if (sumshared/pt->GetNumberOfClusters()>factor){
2104 for (Int_t j=0;j<4;j++){
2105 if (pt->fOverlapLabels[3*j]==0) continue;
2106 if (pt->fOverlapLabels[3*j+1]<5) continue;
2107 if (pt->fRemoval==removalindex) continue;
2108 AliTPCseed * pt2 = (AliTPCseed*)arr->UncheckedAt(pt->fOverlapLabels[3*j+2]);
2110 if (pt2->GetSigma2C()<pt->GetSigma2C()){
2111 // pt->fRemoval = removalindex;
2112 delete arr->RemoveAt(i);
2120 Info("RemoveDouble","Number of tracks after double removal- %d\n",arr->GetEntries());
2129 void AliTPCtrackerMI::SortTracks(TObjArray * arr, Int_t mode) const
2132 //sort tracks in array according mode criteria
2133 Int_t nseed = arr->GetEntriesFast();
2134 for (Int_t i=0; i<nseed; i++) {
2135 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2145 void AliTPCtrackerMI::RemoveUsed(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t removalindex)
2148 //Loop over all tracks and remove "overlaps"
2151 Int_t nseed = arr->GetEntriesFast();
2154 for (Int_t i=0; i<nseed; i++) {
2155 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2157 delete arr->RemoveAt(i);
2161 pt->fBSigned = kFALSE;
2165 nseed = arr->GetEntriesFast();
2172 for (Int_t i=0; i<nseed; i++) {
2173 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2177 Int_t found,foundable,shared;
2179 pt->GetClusterStatistic(0,160,found, foundable,shared,kFALSE);
2181 pt->GetClusterStatistic(0,160,found, foundable,shared,kTRUE);
2183 Double_t factor = factor2;
2184 if (pt->fBConstrain) factor = factor1;
2186 if ((Float_t(shared)/Float_t(found))>factor){
2187 pt->Desactivate(removalindex);
2192 for (Int_t i=0; i<160; i++) {
2193 Int_t index=pt->GetClusterIndex2(i);
2194 if (index<0 || index&0x8000 ) continue;
2195 AliTPCclusterMI *c= pt->fClusterPointer[i];
2197 // if (!c->IsUsed(10)) c->Use(10);
2198 //if (pt->IsActive())
2207 Info("RemoveUsed","\n*****\nNumber of good tracks after shared removal\t%d\n",fNtracks);
2212 void AliTPCtrackerMI::RemoveUsed2(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t minimal)
2215 //Loop over all tracks and remove "overlaps"
2220 Int_t nseed = arr->GetEntriesFast();
2221 Float_t * quality = new Float_t[nseed];
2222 Int_t * indexes = new Int_t[nseed];
2226 for (Int_t i=0; i<nseed; i++) {
2227 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2232 pt->UpdatePoints(); //select first last max dens points
2233 Float_t * points = pt->GetPoints();
2234 if (points[3]<0.8) quality[i] =-1;
2236 quality[i] = (points[2]-points[0])+pt->GetNumberOfClusters();
2238 TMath::Sort(nseed,quality,indexes);
2241 for (Int_t itrack=0; itrack<nseed; itrack++) {
2242 Int_t trackindex = indexes[itrack];
2243 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(trackindex);
2244 if (quality[trackindex]<0){
2246 delete arr->RemoveAt(trackindex);
2249 arr->RemoveAt(trackindex);
2254 Int_t first = Int_t(pt->GetPoints()[0]);
2255 Int_t last = Int_t(pt->GetPoints()[2]);
2256 Double_t factor = (pt->fBConstrain) ? factor1: factor2;
2258 Int_t found,foundable,shared;
2259 pt->GetClusterStatistic(first,last, found, foundable,shared,kFALSE);
2260 Float_t sharedfactor = Float_t(shared+1)/Float_t(found+1);
2261 Bool_t itsgold =kFALSE;
2264 if (pt->fEsd->GetITSclusters(dummy)>4) itsgold= kTRUE;
2268 if (Float_t(shared+1)/Float_t(found+1)>factor){
2269 if (pt->GetKinkIndexes()[0]!=0) continue; //don't remove tracks - part of the kinks
2270 delete arr->RemoveAt(trackindex);
2273 if (pt->GetNumberOfClusters()<50&&(found-0.5*shared)<minimal){ //remove short tracks
2274 if (pt->GetKinkIndexes()[0]!=0) continue; //don't remove tracks - part of the kinks
2275 delete arr->RemoveAt(trackindex);
2281 if (sharedfactor>0.4) continue;
2282 if (pt->GetKinkIndexes()[0]>0) continue;
2283 for (Int_t i=first; i<last; i++) {
2284 Int_t index=pt->GetClusterIndex2(i);
2285 // if (index<0 || index&0x8000 ) continue;
2286 if (index<0 || index&0x8000 ) continue;
2287 AliTPCclusterMI *c= pt->fClusterPointer[i];
2294 Info("RemoveUsed","\n*****\nNumber of good tracks after shared removal\t%d\n",fNtracks);
2300 void AliTPCtrackerMI::UnsignClusters()
2303 // loop over all clusters and unsign them
2306 for (Int_t sec=0;sec<fkNIS;sec++){
2307 for (Int_t row=0;row<fInnerSec->GetNRows();row++){
2308 AliTPCclusterMI *cl = fInnerSec[sec][row].fClusters1;
2309 for (Int_t icl =0;icl< fInnerSec[sec][row].fN1;icl++)
2310 // if (cl[icl].IsUsed(10))
2312 cl = fInnerSec[sec][row].fClusters2;
2313 for (Int_t icl =0;icl< fInnerSec[sec][row].fN2;icl++)
2314 //if (cl[icl].IsUsed(10))
2319 for (Int_t sec=0;sec<fkNOS;sec++){
2320 for (Int_t row=0;row<fOuterSec->GetNRows();row++){
2321 AliTPCclusterMI *cl = fOuterSec[sec][row].fClusters1;
2322 for (Int_t icl =0;icl< fOuterSec[sec][row].fN1;icl++)
2323 //if (cl[icl].IsUsed(10))
2325 cl = fOuterSec[sec][row].fClusters2;
2326 for (Int_t icl =0;icl< fOuterSec[sec][row].fN2;icl++)
2327 //if (cl[icl].IsUsed(10))
2336 void AliTPCtrackerMI::SignClusters(TObjArray * arr, Float_t fnumber, Float_t fdensity)
2339 //sign clusters to be "used"
2341 // snumber and sdensity sign number of sigmas - bellow mean value to be accepted
2342 // loop over "primaries"
2356 Int_t nseed = arr->GetEntriesFast();
2357 for (Int_t i=0; i<nseed; i++) {
2358 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2362 if (!(pt->IsActive())) continue;
2363 Float_t dens = pt->GetNumberOfClusters()/Float_t(pt->fNFoundable);
2364 if ( (dens>0.7) && (pt->GetNumberOfClusters()>70)){
2366 sumdens2+= dens*dens;
2367 sumn += pt->GetNumberOfClusters();
2368 sumn2 += pt->GetNumberOfClusters()*pt->GetNumberOfClusters();
2369 Float_t chi2 = pt->GetChi2()/pt->GetNumberOfClusters();
2372 sumchi2 +=chi2*chi2;
2377 Float_t mdensity = 0.9;
2378 Float_t meann = 130;
2379 Float_t meanchi = 1;
2380 Float_t sdensity = 0.1;
2381 Float_t smeann = 10;
2382 Float_t smeanchi =0.4;
2386 mdensity = sumdens/sum;
2388 meanchi = sumchi/sum;
2390 sdensity = sumdens2/sum-mdensity*mdensity;
2391 sdensity = TMath::Sqrt(sdensity);
2393 smeann = sumn2/sum-meann*meann;
2394 smeann = TMath::Sqrt(smeann);
2396 smeanchi = sumchi2/sum - meanchi*meanchi;
2397 smeanchi = TMath::Sqrt(smeanchi);
2401 //REMOVE SHORT DELTAS or tracks going out of sensitive volume of TPC
2403 for (Int_t i=0; i<nseed; i++) {
2404 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2408 if (pt->fBSigned) continue;
2409 if (pt->fBConstrain) continue;
2410 //if (!(pt->IsActive())) continue;
2412 Int_t found,foundable,shared;
2413 pt->GetClusterStatistic(0,160,found, foundable,shared);
2414 if (shared/float(found)>0.3) {
2415 if (shared/float(found)>0.9 ){
2416 //delete arr->RemoveAt(i);
2421 Bool_t isok =kFALSE;
2422 if ( (pt->fNShared/pt->GetNumberOfClusters()<0.5) &&pt->GetNumberOfClusters()>60)
2424 if ((TMath::Abs(1/pt->GetC())<100.) && (pt->fNShared/pt->GetNumberOfClusters()<0.7))
2426 if (TMath::Abs(pt->GetZ()/pt->GetX())>1.1)
2428 if ( (TMath::Abs(pt->GetSnp()>0.7) && pt->GetD(0,0)>60.))
2432 for (Int_t i=0; i<160; i++) {
2433 Int_t index=pt->GetClusterIndex2(i);
2434 if (index<0) continue;
2435 AliTPCclusterMI *c= pt->fClusterPointer[i];
2437 //if (!(c->IsUsed(10))) c->Use();
2444 Double_t maxchi = meanchi+2.*smeanchi;
2446 for (Int_t i=0; i<nseed; i++) {
2447 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2451 //if (!(pt->IsActive())) continue;
2452 if (pt->fBSigned) continue;
2453 Double_t chi = pt->GetChi2()/pt->GetNumberOfClusters();
2454 if (chi>maxchi) continue;
2457 Float_t dens = pt->GetNumberOfClusters()/Float_t(pt->fNFoundable);
2459 //sign only tracks with enoug big density at the beginning
2461 if ((pt->GetDensityFirst(40)<0.75) && pt->GetNumberOfClusters()<meann) continue;
2464 Double_t mindens = TMath::Max(double(mdensity-sdensity*fdensity*bfactor),0.65);
2465 Double_t minn = TMath::Max(Int_t(meann-fnumber*smeann*bfactor),50);
2467 // if (pt->fBConstrain) mindens = TMath::Max(mdensity-sdensity*fdensity*bfactor,0.65);
2468 if ( (pt->fRemoval==10) && (pt->GetSnp()>0.8)&&(dens>mindens))
2471 if ((dens>mindens && pt->GetNumberOfClusters()>minn) && chi<maxchi ){
2472 //Int_t noc=pt->GetNumberOfClusters();
2473 pt->fBSigned = kTRUE;
2474 for (Int_t i=0; i<160; i++) {
2476 Int_t index=pt->GetClusterIndex2(i);
2477 if (index<0) continue;
2478 AliTPCclusterMI *c= pt->fClusterPointer[i];
2480 // if (!(c->IsUsed(10))) c->Use();
2485 // gLastCheck = nseed;
2493 void AliTPCtrackerMI::StopNotActive(TObjArray * arr, Int_t row0, Float_t th0, Float_t th1, Float_t th2) const
2495 // stop not active tracks
2496 // take th1 as threshold for number of founded to number of foundable on last 10 active rows
2497 // take th2 as threshold for number of founded to number of foundable on last 20 active rows
2498 Int_t nseed = arr->GetEntriesFast();
2500 for (Int_t i=0; i<nseed; i++) {
2501 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2505 if (!(pt->IsActive())) continue;
2506 StopNotActive(pt,row0,th0, th1,th2);
2512 void AliTPCtrackerMI::StopNotActive(AliTPCseed * seed, Int_t row0, Float_t th0, Float_t th1,
2515 // stop not active tracks
2516 // take th1 as threshold for number of founded to number of foundable on last 10 active rows
2517 // take th2 as threshold for number of founded to number of foundable on last 20 active rows
2520 Int_t foundable = 0;
2521 Int_t maxindex = seed->fLastPoint; //last foundable row
2522 if (seed->fNFoundable*th0 > seed->GetNumberOfClusters()) {
2523 seed->Desactivate(10) ;
2527 for (Int_t i=row0; i<maxindex; i++){
2528 Int_t index = seed->GetClusterIndex2(i);
2529 if (index!=-1) foundable++;
2531 if (foundable<=30) sumgood1++;
2532 if (foundable<=50) {
2539 if (foundable>=30.){
2540 if (sumgood1<(th1*30.)) seed->Desactivate(10);
2543 if (sumgood2<(th2*50.)) seed->Desactivate(10);
2547 Int_t AliTPCtrackerMI::RefitInward(AliESD *event)
2550 // back propagation of ESD tracks
2556 //PrepareForProlongation(fSeeds,1);
2557 PropagateForward2(fSeeds);
2560 Int_t nseed = fSeeds->GetEntriesFast();
2561 for (Int_t i=0;i<nseed;i++){
2562 AliTPCseed * seed = (AliTPCseed*) fSeeds->UncheckedAt(i);
2563 if (!seed) continue;
2564 if (seed->GetKinkIndex(0)>0) UpdateKinkQualityD(seed); // update quality informations for kinks
2566 seed->PropagateTo(fParam->GetInnerRadiusLow());
2567 seed->UpdatePoints();
2568 AliESDtrack *esd=event->GetTrack(i);
2569 seed->CookdEdx(0.02,0.6);
2570 CookLabel(seed,0.1); //For comparison only
2572 if (AliTPCReconstructor::StreamLevel()>0 && seed!=0&&esd!=0) {
2573 TTreeSRedirector &cstream = *fDebugStreamer;
2579 if (seed->GetNumberOfClusters()>15){
2580 esd->UpdateTrackParams(seed,AliESDtrack::kTPCrefit);
2581 esd->SetTPCPoints(seed->GetPoints());
2582 esd->SetTPCPointsF(seed->fNFoundable);
2583 Int_t ndedx = seed->fNCDEDX[0]+seed->fNCDEDX[1]+seed->fNCDEDX[2]+seed->fNCDEDX[3];
2584 Float_t sdedx = (seed->fSDEDX[0]+seed->fSDEDX[1]+seed->fSDEDX[2]+seed->fSDEDX[3])*0.25;
2585 Float_t dedx = seed->GetdEdx();
2586 esd->SetTPCsignal(dedx, sdedx, ndedx);
2590 //printf("problem\n");
2593 //FindKinks(fSeeds,event);
2594 Info("RefitInward","Number of refitted tracks %d",ntracks);
2601 Int_t AliTPCtrackerMI::PropagateBack(AliESD *event)
2604 // back propagation of ESD tracks
2610 PropagateBack(fSeeds);
2611 RemoveUsed2(fSeeds,0.4,0.4,20);
2613 Int_t nseed = fSeeds->GetEntriesFast();
2615 for (Int_t i=0;i<nseed;i++){
2616 AliTPCseed * seed = (AliTPCseed*) fSeeds->UncheckedAt(i);
2617 if (!seed) continue;
2618 if (seed->GetKinkIndex(0)<0) UpdateKinkQualityM(seed); // update quality informations for kinks
2619 seed->UpdatePoints();
2620 AliESDtrack *esd=event->GetTrack(i);
2621 seed->CookdEdx(0.02,0.6);
2622 CookLabel(seed,0.1); //For comparison only
2623 if (seed->GetNumberOfClusters()>15){
2624 esd->UpdateTrackParams(seed,AliESDtrack::kTPCout);
2625 esd->SetTPCPoints(seed->GetPoints());
2626 esd->SetTPCPointsF(seed->fNFoundable);
2627 Int_t ndedx = seed->fNCDEDX[0]+seed->fNCDEDX[1]+seed->fNCDEDX[2]+seed->fNCDEDX[3];
2628 Float_t sdedx = (seed->fSDEDX[0]+seed->fSDEDX[1]+seed->fSDEDX[2]+seed->fSDEDX[3])*0.25;
2629 Float_t dedx = seed->GetdEdx();
2630 esd->SetTPCsignal(dedx, sdedx, ndedx);
2632 Int_t eventnumber = event->GetEventNumber();// patch 28 fev 06
2633 (*fDebugStreamer)<<"Cback"<<
2635 "EventNr="<<eventnumber<<
2636 "\n"; // patch 28 fev 06
2639 //FindKinks(fSeeds,event);
2640 Info("PropagateBack","Number of back propagated tracks %d",ntracks);
2647 void AliTPCtrackerMI::DeleteSeeds()
2651 Int_t nseed = fSeeds->GetEntriesFast();
2652 for (Int_t i=0;i<nseed;i++){
2653 AliTPCseed * seed = (AliTPCseed*)fSeeds->At(i);
2654 if (seed) delete fSeeds->RemoveAt(i);
2660 void AliTPCtrackerMI::ReadSeeds(AliESD *event, Int_t direction)
2663 //read seeds from the event
2665 Int_t nentr=event->GetNumberOfTracks();
2667 Info("ReadSeeds", "Number of ESD tracks: %d\n", nentr);
2672 fSeeds = new TObjArray(nentr);
2676 for (Int_t i=0; i<nentr; i++) {
2677 AliESDtrack *esd=event->GetTrack(i);
2678 ULong_t status=esd->GetStatus();
2679 if (!(status&AliESDtrack::kTPCin)) continue;
2680 AliTPCtrack t(*esd);
2681 t.SetNumberOfClusters(0);
2682 // AliTPCseed *seed = new AliTPCseed(t,t.GetAlpha());
2683 AliTPCseed *seed = new AliTPCseed(t/*,t.GetAlpha()*/);
2684 for (Int_t ikink=0;ikink<3;ikink++) {
2685 Int_t index = esd->GetKinkIndex(ikink);
2686 seed->GetKinkIndexes()[ikink] = index;
2687 if (index==0) continue;
2688 index = TMath::Abs(index);
2689 AliESDkink * kink = fEvent->GetKink(index-1);
2690 if (kink&&esd->GetKinkIndex(ikink)<0){
2691 if ((status & AliESDtrack::kTRDrefit) != 0) kink->SetStatus(1,2);
2692 if ((status & AliESDtrack::kITSout) != 0) kink->SetStatus(1,0);
2694 if (kink&&esd->GetKinkIndex(ikink)>0){
2695 if ((status & AliESDtrack::kTRDrefit) != 0) kink->SetStatus(1,6);
2696 if ((status & AliESDtrack::kITSout) != 0) kink->SetStatus(1,4);
2700 if (((status&AliESDtrack::kITSout)==0)&&(direction==1)) seed->ResetCovariance();
2701 if ( direction ==2 &&(status & AliESDtrack::kTRDrefit) == 0 ) seed->ResetCovariance();
2702 if ( direction ==2 && ((status & AliESDtrack::kTPCout) == 0) ) {
2707 if ( direction ==2 &&(status & AliESDtrack::kTRDrefit) > 0 ) {
2708 Double_t par0[5],par1[5],alpha,x;
2709 esd->GetInnerExternalParameters(alpha,x,par0);
2710 esd->GetExternalParameters(x,par1);
2711 Double_t delta1 = TMath::Abs(par0[4]-par1[4])/(0.000000001+TMath::Abs(par0[4]+par1[4]));
2712 Double_t delta2 = TMath::Abs(par0[3]-par1[3]);
2714 if (esd->GetTRDncls()>0) trdchi2 = esd->GetTRDchi2()/esd->GetTRDncls();
2715 //reset covariance if suspicious
2716 if ( (delta1>0.1) || (delta2>0.006) ||trdchi2>7.)
2717 seed->ResetCovariance();
2722 // rotate to the local coordinate system
2724 fSectors=fInnerSec; fN=fkNIS;
2725 Double_t alpha=seed->GetAlpha() - fSectors->GetAlphaShift();
2726 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
2727 if (alpha < 0. ) alpha += 2.*TMath::Pi();
2728 Int_t ns=Int_t(alpha/fSectors->GetAlpha())%fN;
2729 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
2730 if (alpha<-TMath::Pi()) alpha += 2*TMath::Pi();
2731 if (alpha>=TMath::Pi()) alpha -= 2*TMath::Pi();
2732 alpha-=seed->GetAlpha();
2733 if (!seed->Rotate(alpha)) {
2739 if (esd->GetKinkIndex(0)<=0){
2740 for (Int_t irow=0;irow<160;irow++){
2741 Int_t index = seed->GetClusterIndex2(irow);
2744 AliTPCclusterMI * cl = GetClusterMI(index);
2745 seed->fClusterPointer[irow] = cl;
2747 if ((index & 0x8000)==0){
2748 cl->Use(10); // accepted cluster
2750 cl->Use(6); // close cluster not accepted
2753 Info("ReadSeeds","Not found cluster");
2758 fSeeds->AddAt(seed,i);
2764 //_____________________________________________________________________________
2765 void AliTPCtrackerMI::MakeSeeds3(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4],
2766 Float_t deltay, Int_t ddsec) {
2767 //-----------------------------------------------------------------
2768 // This function creates track seeds.
2769 // SEEDING WITH VERTEX CONSTRAIN
2770 //-----------------------------------------------------------------
2771 // cuts[0] - fP4 cut
2772 // cuts[1] - tan(phi) cut
2773 // cuts[2] - zvertex cut
2774 // cuts[3] - fP3 cut
2782 Double_t x[5], c[15];
2783 // Int_t di = i1-i2;
2785 AliTPCseed * seed = new AliTPCseed;
2786 Double_t alpha=fSectors->GetAlpha(), shift=fSectors->GetAlphaShift();
2787 Double_t cs=cos(alpha), sn=sin(alpha);
2789 // Double_t x1 =fOuterSec->GetX(i1);
2790 //Double_t xx2=fOuterSec->GetX(i2);
2792 Double_t x1 =GetXrow(i1);
2793 Double_t xx2=GetXrow(i2);
2795 Double_t x3=GetX(), y3=GetY(), z3=GetZ();
2797 Int_t imiddle = (i2+i1)/2; //middle pad row index
2798 Double_t xm = GetXrow(imiddle); // radius of middle pad-row
2799 const AliTPCRow& krm=GetRow(sec,imiddle); //middle pad -row
2803 const AliTPCRow& kr1=GetRow(ns,i1);
2804 Double_t ymax = GetMaxY(i1)-kr1.fDeadZone-1.5;
2805 Double_t ymaxm = GetMaxY(imiddle)-kr1.fDeadZone-1.5;
2808 // change cut on curvature if it can't reach this layer
2809 // maximal curvature set to reach it
2810 Double_t dvertexmax = TMath::Sqrt((x1-x3)*(x1-x3)+(ymax+5-y3)*(ymax+5-y3));
2811 if (dvertexmax*0.5*cuts[0]>0.85){
2812 cuts[0] = 0.85/(dvertexmax*0.5+1.);
2814 Double_t r2min = 1/(cuts[0]*cuts[0]); //minimal square of radius given by cut
2817 if (deltay>0) ddsec = 0;
2818 // loop over clusters
2819 for (Int_t is=0; is < kr1; is++) {
2821 if (kr1[is]->IsUsed(10)) continue;
2822 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
2823 //if (TMath::Abs(y1)>ymax) continue;
2825 if (deltay>0 && TMath::Abs(ymax-TMath::Abs(y1))> deltay ) continue; // seed only at the edge
2827 // find possible directions
2828 Float_t anglez = (z1-z3)/(x1-x3);
2829 Float_t extraz = z1 - anglez*(x1-xx2); // extrapolated z
2832 //find rotation angles relative to line given by vertex and point 1
2833 Double_t dvertex2 = (x1-x3)*(x1-x3)+(y1-y3)*(y1-y3);
2834 Double_t dvertex = TMath::Sqrt(dvertex2);
2835 Double_t angle13 = TMath::ATan((y1-y3)/(x1-x3));
2836 Double_t cs13 = cos(-angle13), sn13 = sin(-angle13);
2839 // loop over 2 sectors
2845 Double_t dddz1=0; // direction of delta inclination in z axis
2852 for (Int_t dsec = dsec1; dsec<=dsec2;dsec++){
2853 Int_t sec2 = sec + dsec;
2855 // AliTPCRow& kr2 = fOuterSec[(sec2+fkNOS)%fkNOS][i2];
2856 //AliTPCRow& kr2m = fOuterSec[(sec2+fkNOS)%fkNOS][imiddle];
2857 AliTPCRow& kr2 = GetRow((sec2+fkNOS)%fkNOS,i2);
2858 AliTPCRow& kr2m = GetRow((sec2+fkNOS)%fkNOS,imiddle);
2859 Int_t index1 = TMath::Max(kr2.Find(extraz-0.6-dddz1*TMath::Abs(z1)*0.05)-1,0);
2860 Int_t index2 = TMath::Min(kr2.Find(extraz+0.6+dddz2*TMath::Abs(z1)*0.05)+1,kr2);
2862 // rotation angles to p1-p3
2863 Double_t cs13r = cos(-angle13+dsec*alpha)/dvertex, sn13r = sin(-angle13+dsec*alpha)/dvertex;
2864 Double_t x2, y2, z2;
2866 // Double_t dymax = maxangle*TMath::Abs(x1-xx2);
2869 Double_t dxx0 = (xx2-x3)*cs13r;
2870 Double_t dyy0 = (xx2-x3)*sn13r;
2871 for (Int_t js=index1; js < index2; js++) {
2872 const AliTPCclusterMI *kcl = kr2[js];
2873 if (kcl->IsUsed(10)) continue;
2875 //calcutate parameters
2877 Double_t yy0 = dyy0 +(kcl->GetY()-y3)*cs13r;
2879 if (TMath::Abs(yy0)<0.000001) continue;
2880 Double_t xx0 = dxx0 -(kcl->GetY()-y3)*sn13r;
2881 Double_t y0 = 0.5*(xx0*xx0+yy0*yy0-xx0)/yy0;
2882 Double_t r02 = (0.25+y0*y0)*dvertex2;
2883 //curvature (radius) cut
2884 if (r02<r2min) continue;
2888 Double_t c0 = 1/TMath::Sqrt(r02);
2892 //Double_t dfi0 = 2.*TMath::ASin(dvertex*c0*0.5);
2893 //Double_t dfi1 = 2.*TMath::ASin(TMath::Sqrt(yy0*yy0+(1-xx0)*(1-xx0))*dvertex*c0*0.5);
2894 Double_t dfi0 = 2.*AliTPCFastMath::FastAsin(dvertex*c0*0.5);
2895 Double_t dfi1 = 2.*AliTPCFastMath::FastAsin(TMath::Sqrt(yy0*yy0+(1-xx0)*(1-xx0))*dvertex*c0*0.5);
2898 Double_t z0 = kcl->GetZ();
2899 Double_t zzzz2 = z1-(z1-z3)*dfi1/dfi0;
2900 if (TMath::Abs(zzzz2-z0)>0.5) continue;
2903 Double_t dip = (z1-z0)*c0/dfi1;
2904 Double_t x0 = (0.5*cs13+y0*sn13)*dvertex*c0;
2915 x2= xx2*cs-y2*sn*dsec;
2916 y2=+xx2*sn*dsec+y2*cs;
2926 // do we have cluster at the middle ?
2928 GetProlongation(x1,xm,x,ym,zm);
2930 AliTPCclusterMI * cm=0;
2931 if (TMath::Abs(ym)-ymaxm<0){
2932 cm = krm.FindNearest2(ym,zm,1.0,0.6,dummy);
2933 if ((!cm) || (cm->IsUsed(10))) {
2938 // rotate y1 to system 0
2939 // get state vector in rotated system
2940 Double_t yr1 = (-0.5*sn13+y0*cs13)*dvertex*c0;
2941 Double_t xr2 = x0*cs+yr1*sn*dsec;
2942 Double_t xr[5]={kcl->GetY(),kcl->GetZ(), xr2, dip, c0};
2944 GetProlongation(xx2,xm,xr,ym,zm);
2945 if (TMath::Abs(ym)-ymaxm<0){
2946 cm = kr2m.FindNearest2(ym,zm,1.0,0.6,dummy);
2947 if ((!cm) || (cm->IsUsed(10))) {
2957 dym = ym - cm->GetY();
2958 dzm = zm - cm->GetZ();
2965 Double_t sy1=kr1[is]->GetSigmaY2()*2., sz1=kr1[is]->GetSigmaZ2()*2.;
2966 Double_t sy2=kcl->GetSigmaY2()*2., sz2=kcl->GetSigmaZ2()*2.;
2967 //Double_t sy3=400*3./12., sy=0.1, sz=0.1;
2968 Double_t sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
2969 //Double_t sy3=25000*x[4]*x[4]*60+0.5, sy=0.1, sz=0.1;
2971 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
2972 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
2973 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
2974 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
2975 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
2976 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
2978 Double_t f30=(F3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
2979 Double_t f31=(F3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
2980 Double_t f32=(F3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
2981 Double_t f34=(F3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
2985 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
2986 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
2987 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
2988 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
2989 c[13]=f30*sy1*f40+f32*sy2*f42;
2990 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
2992 // if (!BuildSeed(kr1[is],kcl,0,x1,x2,x3,x,c)) continue;
2994 UInt_t index=kr1.GetIndex(is);
2995 AliTPCseed *track=new(seed) AliTPCseed(index, x, c, x1, ns*alpha+shift);
2997 track->fIsSeeding = kTRUE;
3004 FollowProlongation(*track, (i1+i2)/2,1);
3005 Int_t foundable,found,shared;
3006 track->GetClusterStatistic((i1+i2)/2,i1, found, foundable, shared, kTRUE);
3007 if ((found<0.55*foundable) || shared>0.5*found || (track->GetSigmaY2()+track->GetSigmaZ2())>0.5){
3009 seed->~AliTPCseed();
3015 FollowProlongation(*track, i2,1);
3019 track->fBConstrain =1;
3020 // track->fLastPoint = i1+fInnerSec->GetNRows(); // first cluster in track position
3021 track->fLastPoint = i1; // first cluster in track position
3022 track->fFirstPoint = track->fLastPoint;
3024 if (track->GetNumberOfClusters()<(i1-i2)*0.5 ||
3025 track->GetNumberOfClusters() < track->fNFoundable*0.6 ||
3026 track->fNShared>0.4*track->GetNumberOfClusters() ) {
3028 seed->~AliTPCseed();
3032 // Z VERTEX CONDITION
3034 zv = track->GetZ()+track->GetTgl()/track->GetC()*
3035 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
3036 if (TMath::Abs(zv-z3)>cuts[2]) {
3037 FollowProlongation(*track, TMath::Max(i2-20,0));
3038 zv = track->GetZ()+track->GetTgl()/track->GetC()*
3039 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
3040 if (TMath::Abs(zv-z3)>cuts[2]){
3041 FollowProlongation(*track, TMath::Max(i2-40,0));
3042 zv = track->GetZ()+track->GetTgl()/track->GetC()*
3043 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
3044 if (TMath::Abs(zv-z3)>cuts[2] &&(track->GetNumberOfClusters() > track->fNFoundable*0.7)){
3045 // make seed without constrain
3046 AliTPCseed * track2 = MakeSeed(track,0.2,0.5,1.);
3047 FollowProlongation(*track2, i2,1);
3048 track2->fBConstrain = kFALSE;
3049 track2->fSeedType = 1;
3050 arr->AddLast(track2);
3052 seed->~AliTPCseed();
3057 seed->~AliTPCseed();
3064 track->fSeedType =0;
3065 arr->AddLast(track);
3066 seed = new AliTPCseed;
3068 // don't consider other combinations
3069 if (track->GetNumberOfClusters() > track->fNFoundable*0.8)
3075 Info("MakeSeeds3","\nSeeding statistic:\t%d\t%d\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin,nout1,nout2);
3081 void AliTPCtrackerMI::MakeSeeds5(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4],
3086 //-----------------------------------------------------------------
3087 // This function creates track seeds.
3088 //-----------------------------------------------------------------
3089 // cuts[0] - fP4 cut
3090 // cuts[1] - tan(phi) cut
3091 // cuts[2] - zvertex cut
3092 // cuts[3] - fP3 cut
3102 Double_t x[5], c[15];
3104 // make temporary seed
3105 AliTPCseed * seed = new AliTPCseed;
3106 Double_t alpha=fOuterSec->GetAlpha(), shift=fOuterSec->GetAlphaShift();
3107 // Double_t cs=cos(alpha), sn=sin(alpha);
3112 Double_t x1 = GetXrow(i1-1);
3113 const AliTPCRow& kr1=GetRow(sec,i1-1);
3114 Double_t y1max = GetMaxY(i1-1)-kr1.fDeadZone-1.5;
3116 Double_t x1p = GetXrow(i1);
3117 const AliTPCRow& kr1p=GetRow(sec,i1);
3119 Double_t x1m = GetXrow(i1-2);
3120 const AliTPCRow& kr1m=GetRow(sec,i1-2);
3123 //last 3 padrow for seeding
3124 AliTPCRow& kr3 = GetRow((sec+fkNOS)%fkNOS,i1-7);
3125 Double_t x3 = GetXrow(i1-7);
3126 // Double_t y3max= GetMaxY(i1-7)-kr3.fDeadZone-1.5;
3128 AliTPCRow& kr3p = GetRow((sec+fkNOS)%fkNOS,i1-6);
3129 Double_t x3p = GetXrow(i1-6);
3131 AliTPCRow& kr3m = GetRow((sec+fkNOS)%fkNOS,i1-8);
3132 Double_t x3m = GetXrow(i1-8);
3137 Int_t im = i1-4; //middle pad row index
3138 Double_t xm = GetXrow(im); // radius of middle pad-row
3139 const AliTPCRow& krm=GetRow(sec,im); //middle pad -row
3140 // Double_t ymmax = GetMaxY(im)-kr1.fDeadZone-1.5;
3143 Double_t deltax = x1-x3;
3144 Double_t dymax = deltax*cuts[1];
3145 Double_t dzmax = deltax*cuts[3];
3147 // loop over clusters
3148 for (Int_t is=0; is < kr1; is++) {
3150 if (kr1[is]->IsUsed(10)) continue;
3151 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
3153 if (deltay>0 && TMath::Abs(y1max-TMath::Abs(y1))> deltay ) continue; // seed only at the edge
3155 Int_t index1 = TMath::Max(kr3.Find(z1-dzmax)-1,0);
3156 Int_t index2 = TMath::Min(kr3.Find(z1+dzmax)+1,kr3);
3162 for (Int_t js=index1; js < index2; js++) {
3163 const AliTPCclusterMI *kcl = kr3[js];
3164 if (kcl->IsUsed(10)) continue;
3166 // apply angular cuts
3167 if (TMath::Abs(y1-y3)>dymax) continue;
3170 if (TMath::Abs(z1-z3)>dzmax) continue;
3172 Double_t angley = (y1-y3)/(x1-x3);
3173 Double_t anglez = (z1-z3)/(x1-x3);
3175 Double_t erry = TMath::Abs(angley)*(x1-x1m)*0.5+0.5;
3176 Double_t errz = TMath::Abs(anglez)*(x1-x1m)*0.5+0.5;
3178 Double_t yyym = angley*(xm-x1)+y1;
3179 Double_t zzzm = anglez*(xm-x1)+z1;
3181 const AliTPCclusterMI *kcm = krm.FindNearest2(yyym,zzzm,erry,errz,index);
3183 if (kcm->IsUsed(10)) continue;
3185 erry = TMath::Abs(angley)*(x1-x1m)*0.4+0.5;
3186 errz = TMath::Abs(anglez)*(x1-x1m)*0.4+0.5;
3193 // look around first
3194 const AliTPCclusterMI *kc1m = kr1m.FindNearest2(angley*(x1m-x1)+y1,
3200 if (kc1m->IsUsed(10)) used++;
3202 const AliTPCclusterMI *kc1p = kr1p.FindNearest2(angley*(x1p-x1)+y1,
3208 if (kc1p->IsUsed(10)) used++;
3210 if (used>1) continue;
3211 if (found<1) continue;
3215 const AliTPCclusterMI *kc3m = kr3m.FindNearest2(angley*(x3m-x3)+y3,
3221 if (kc3m->IsUsed(10)) used++;
3225 const AliTPCclusterMI *kc3p = kr3p.FindNearest2(angley*(x3p-x3)+y3,
3231 if (kc3p->IsUsed(10)) used++;
3235 if (used>1) continue;
3236 if (found<3) continue;
3246 x[4]=F1(x1,y1,x2,y2,x3,y3);
3247 //if (TMath::Abs(x[4]) >= cuts[0]) continue;
3250 x[2]=F2(x1,y1,x2,y2,x3,y3);
3253 x[3]=F3n(x1,y1,x2,y2,z1,z2,x[4]);
3254 //if (TMath::Abs(x[3]) > cuts[3]) continue;
3258 Double_t sy1=0.1, sz1=0.1;
3259 Double_t sy2=0.1, sz2=0.1;
3260 Double_t sy3=0.1, sy=0.1, sz=0.1;
3262 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
3263 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
3264 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
3265 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
3266 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
3267 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
3269 Double_t f30=(F3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
3270 Double_t f31=(F3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
3271 Double_t f32=(F3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
3272 Double_t f34=(F3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
3276 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3277 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3278 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3279 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3280 c[13]=f30*sy1*f40+f32*sy2*f42;
3281 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3283 // if (!BuildSeed(kr1[is],kcl,0,x1,x2,x3,x,c)) continue;
3285 UInt_t index=kr1.GetIndex(is);
3286 AliTPCseed *track=new(seed) AliTPCseed(index, x, c, x1, sec*alpha+shift);
3288 track->fIsSeeding = kTRUE;
3291 FollowProlongation(*track, i1-7,1);
3292 if (track->GetNumberOfClusters() < track->fNFoundable*0.75 ||
3293 track->fNShared>0.6*track->GetNumberOfClusters() || ( track->GetSigmaY2()+ track->GetSigmaZ2())>0.6){
3295 seed->~AliTPCseed();
3301 FollowProlongation(*track, i2,1);
3302 track->fBConstrain =0;
3303 track->fLastPoint = i1+fInnerSec->GetNRows(); // first cluster in track position
3304 track->fFirstPoint = track->fLastPoint;
3306 if (track->GetNumberOfClusters()<(i1-i2)*0.5 ||
3307 track->GetNumberOfClusters()<track->fNFoundable*0.7 ||
3308 track->fNShared>2. || track->GetChi2()/track->GetNumberOfClusters()>6 || ( track->GetSigmaY2()+ track->GetSigmaZ2())>0.5 ) {
3310 seed->~AliTPCseed();
3315 FollowProlongation(*track, TMath::Max(i2-10,0),1);
3316 AliTPCseed * track2 = MakeSeed(track,0.2,0.5,0.9);
3317 FollowProlongation(*track2, i2,1);
3318 track2->fBConstrain = kFALSE;
3319 track2->fSeedType = 4;
3320 arr->AddLast(track2);
3322 seed->~AliTPCseed();
3326 //arr->AddLast(track);
3327 //seed = new AliTPCseed;
3333 Info("MakeSeeds5","\nSeeding statiistic:\t%d\t%d\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin,nout1,nout2,nout3);
3339 //_____________________________________________________________________________
3340 void AliTPCtrackerMI::MakeSeeds2(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t */*cuts[4]*/,
3341 Float_t deltay, Bool_t /*bconstrain*/) {
3342 //-----------------------------------------------------------------
3343 // This function creates track seeds - without vertex constraint
3344 //-----------------------------------------------------------------
3345 // cuts[0] - fP4 cut - not applied
3346 // cuts[1] - tan(phi) cut
3347 // cuts[2] - zvertex cut - not applied
3348 // cuts[3] - fP3 cut
3358 Double_t alpha=fOuterSec->GetAlpha(), shift=fOuterSec->GetAlphaShift();
3359 // Double_t cs=cos(alpha), sn=sin(alpha);
3360 Int_t row0 = (i1+i2)/2;
3361 Int_t drow = (i1-i2)/2;
3362 const AliTPCRow& kr0=fSectors[sec][row0];
3365 AliTPCpolyTrack polytrack;
3366 Int_t nclusters=fSectors[sec][row0];
3367 AliTPCseed * seed = new AliTPCseed;
3372 for (Int_t is=0; is < nclusters; is++) { //LOOP over clusters
3374 Int_t nfoundable =0;
3375 for (Int_t iter =1; iter<2; iter++){ //iterations
3376 const AliTPCRow& krm=fSectors[sec][row0-iter];
3377 const AliTPCRow& krp=fSectors[sec][row0+iter];
3378 const AliTPCclusterMI * cl= kr0[is];
3380 if (cl->IsUsed(10)) {
3386 Double_t x = kr0.GetX();
3387 // Initialization of the polytrack
3392 Double_t y0= cl->GetY();
3393 Double_t z0= cl->GetZ();
3397 Double_t ymax = fSectors->GetMaxY(row0)-kr0.fDeadZone-1.5;
3398 if (deltay>0 && TMath::Abs(ymax-TMath::Abs(y0))> deltay ) continue; // seed only at the edge
3400 erry = (0.5)*cl->GetSigmaY2()/TMath::Sqrt(cl->GetQ())*6;
3401 errz = (0.5)*cl->GetSigmaZ2()/TMath::Sqrt(cl->GetQ())*6;
3402 polytrack.AddPoint(x,y0,z0,erry, errz);
3405 if (cl->IsUsed(10)) sumused++;
3408 Float_t roady = (5*TMath::Sqrt(cl->GetSigmaY2()+0.2)+1.)*iter;
3409 Float_t roadz = (5*TMath::Sqrt(cl->GetSigmaZ2()+0.2)+1.)*iter;
3412 AliTPCclusterMI * cl1 = krm.FindNearest(y0,z0,roady,roadz);
3413 if (cl1 && TMath::Abs(ymax-TMath::Abs(y0))) {
3414 erry = (0.5)*cl1->GetSigmaY2()/TMath::Sqrt(cl1->GetQ())*3;
3415 errz = (0.5)*cl1->GetSigmaZ2()/TMath::Sqrt(cl1->GetQ())*3;
3416 if (cl1->IsUsed(10)) sumused++;
3417 polytrack.AddPoint(x,cl1->GetY(),cl1->GetZ(),erry,errz);
3421 AliTPCclusterMI * cl2 = krp.FindNearest(y0,z0,roady,roadz);
3423 erry = (0.5)*cl2->GetSigmaY2()/TMath::Sqrt(cl2->GetQ())*3;
3424 errz = (0.5)*cl2->GetSigmaZ2()/TMath::Sqrt(cl2->GetQ())*3;
3425 if (cl2->IsUsed(10)) sumused++;
3426 polytrack.AddPoint(x,cl2->GetY(),cl2->GetZ(),erry,errz);
3429 if (sumused>0) continue;
3431 polytrack.UpdateParameters();
3437 nfoundable = polytrack.GetN();
3438 nfound = nfoundable;
3440 for (Int_t ddrow = iter+1; ddrow<drow;ddrow++){
3441 Float_t maxdist = 0.8*(1.+3./(ddrow));
3442 for (Int_t delta = -1;delta<=1;delta+=2){
3443 Int_t row = row0+ddrow*delta;
3444 kr = &(fSectors[sec][row]);
3445 Double_t xn = kr->GetX();
3446 Double_t ymax = fSectors->GetMaxY(row)-kr->fDeadZone-1.5;
3447 polytrack.GetFitPoint(xn,yn,zn);
3448 if (TMath::Abs(yn)>ymax) continue;
3450 AliTPCclusterMI * cln = kr->FindNearest(yn,zn,roady,roadz);
3452 Float_t dist = TMath::Sqrt( (yn-cln->GetY())*(yn-cln->GetY())+(zn-cln->GetZ())*(zn-cln->GetZ()));
3455 erry = (dist+0.3)*cln->GetSigmaY2()/TMath::Sqrt(cln->GetQ())*(1.+1./(ddrow));
3456 errz = (dist+0.3)*cln->GetSigmaZ2()/TMath::Sqrt(cln->GetQ())*(1.+1./(ddrow));
3457 if (cln->IsUsed(10)) {
3458 // printf("used\n");
3466 polytrack.AddPoint(xn,cln->GetY(),cln->GetZ(),erry, errz);
3471 if ( (sumused>3) || (sumused>0.5*nfound) || (nfound<0.6*nfoundable)) break;
3472 polytrack.UpdateParameters();
3475 if ( (sumused>3) || (sumused>0.5*nfound)) {
3476 //printf("sumused %d\n",sumused);
3481 polytrack.GetFitDerivation(kr0.GetX(),dy,dz);
3482 AliTPCpolyTrack track2;
3484 polytrack.Refit(track2,0.5+TMath::Abs(dy)*0.3,0.4+TMath::Abs(dz)*0.3);
3485 if (track2.GetN()<0.5*nfoundable) continue;
3488 if ((nfound>0.6*nfoundable) &&( nfoundable>0.4*(i1-i2))) {
3490 // test seed with and without constrain
3491 for (Int_t constrain=0; constrain<=0;constrain++){
3492 // add polytrack candidate
3494 Double_t x[5], c[15];
3495 Double_t x1,x2,x3,y1,y2,y3,z1,z2,z3;
3496 track2.GetBoundaries(x3,x1);
3498 track2.GetFitPoint(x1,y1,z1);
3499 track2.GetFitPoint(x2,y2,z2);
3500 track2.GetFitPoint(x3,y3,z3);
3502 //is track pointing to the vertex ?
3505 polytrack.GetFitPoint(x0,y0,z0);
3518 x[4]=F1(x1,y1,x2,y2,x3,y3);
3520 // if (TMath::Abs(x[4]) >= cuts[0]) continue; //
3521 x[2]=F2(x1,y1,x2,y2,x3,y3);
3523 //if (TMath::Abs(x[4]*x1-x[2]) >= cuts[1]) continue;
3524 //x[3]=F3(x1,y1,x2,y2,z1,z2);
3525 x[3]=F3n(x1,y1,x3,y3,z1,z3,x[4]);
3526 //if (TMath::Abs(x[3]) > cuts[3]) continue;
3529 Double_t sy =0.1, sz =0.1;
3530 Double_t sy1=0.02, sz1=0.02;
3531 Double_t sy2=0.02, sz2=0.02;
3535 sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
3538 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
3539 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
3540 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
3541 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
3542 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
3543 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
3545 Double_t f30=(F3(x1,y1+sy,x3,y3,z1,z3)-x[3])/sy;
3546 Double_t f31=(F3(x1,y1,x3,y3,z1+sz,z3)-x[3])/sz;
3547 Double_t f32=(F3(x1,y1,x3,y3+sy,z1,z3)-x[3])/sy;
3548 Double_t f34=(F3(x1,y1,x3,y3,z1,z3+sz)-x[3])/sz;
3553 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3554 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3555 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3556 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3557 c[13]=f30*sy1*f40+f32*sy2*f42;
3558 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3560 //Int_t row1 = fSectors->GetRowNumber(x1);
3561 Int_t row1 = GetRowNumber(x1);
3565 AliTPCseed *track=new (seed) AliTPCseed(index, x, c, x1, sec*alpha+shift);
3566 track->fIsSeeding = kTRUE;
3567 Int_t rc=FollowProlongation(*track, i2);
3568 if (constrain) track->fBConstrain =1;
3570 track->fBConstrain =0;
3571 track->fLastPoint = row1+fInnerSec->GetNRows(); // first cluster in track position
3572 track->fFirstPoint = track->fLastPoint;
3574 if (rc==0 || track->GetNumberOfClusters()<(i1-i2)*0.5 ||
3575 track->GetNumberOfClusters() < track->fNFoundable*0.6 ||
3576 track->fNShared>0.4*track->GetNumberOfClusters()) {
3579 seed->~AliTPCseed();
3582 arr->AddLast(track);
3583 seed = new AliTPCseed;
3587 } // if accepted seed
3590 Info("MakeSeeds2","\nSeeding statiistic:\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin3);
3596 AliTPCseed *AliTPCtrackerMI::MakeSeed(AliTPCseed *track, Float_t r0, Float_t r1, Float_t r2)
3600 //reseed using track points
3601 Int_t p0 = int(r0*track->GetNumberOfClusters()); // point 0
3602 Int_t p1 = int(r1*track->GetNumberOfClusters());
3603 Int_t p2 = int(r2*track->GetNumberOfClusters()); // last point
3605 Double_t x0[3],x1[3],x2[3];
3606 for (Int_t i=0;i<3;i++){
3612 // find track position at given ratio of the length
3613 Int_t sec0=0, sec1=0, sec2=0;
3616 for (Int_t i=0;i<160;i++){
3617 if (track->fClusterPointer[i]){