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 -
24 // run AliTPCFindClusters.C macro - clusters neccessary for tracker are founded
25 // run AliTPCFindTracksMI.C macro - to find tracks
26 // tracks are written to AliTPCtracks.root file
27 // for comparison also seeds are written to the same file - to special branch
28 //-------------------------------------------------------
34 #include <TObjArray.h>
37 #include <TClonesArray.h>
39 #include "Riostream.h"
41 #include "AliTPCclusterMI.h"
42 #include "AliComplexCluster.h"
43 #include "AliTPCParam.h"
44 #include "AliTPCClustersRow.h"
45 #include "AliTPCpolyTrack.h"
46 #include "TStopwatch.h"
50 #include "AliRunLoader.h"
52 #include "AliTPCreco.h"
53 #include "AliTPCtrackerMI.h"
58 ClassImp(AliTPCtrackerMI)
61 class AliTPCFastMath {
64 static Double_t FastAsin(Double_t x);
66 static Double_t fgFastAsin[20000]; //lookup table for fast asin computation
69 Double_t AliTPCFastMath::fgFastAsin[20000];
70 AliTPCFastMath gAliTPCFastMath;
72 AliTPCFastMath::AliTPCFastMath(){
74 // initialized lookup table;
75 for (Int_t i=0;i<10000;i++){
76 fgFastAsin[2*i] = TMath::ASin(i/10000.);
77 fgFastAsin[2*i+1] = (TMath::ASin((i+1)/10000.)-fgFastAsin[2*i]);
81 Double_t AliTPCFastMath::FastAsin(Double_t x){
83 // return asin using lookup table
85 Int_t index = int(x*10000);
86 return fgFastAsin[2*index]+(x*10000.-index)*fgFastAsin[2*index+1];
89 Int_t index = int(x*10000);
90 return -(fgFastAsin[2*index]+(x*10000.-index)*fgFastAsin[2*index+1]);
96 Int_t AliTPCtrackerMI::UpdateTrack(AliTPCseed * track, Int_t accept){
98 //update track information using current cluster - track->fCurrentCluster
101 AliTPCclusterMI* c =track->fCurrentCluster;
102 if (accept>0) track->fCurrentClusterIndex1 |=0x8000; //sign not accepted clusters
104 UInt_t i = track->fCurrentClusterIndex1;
106 Int_t sec=(i&0xff000000)>>24;
107 //Int_t row = (i&0x00ff0000)>>16;
108 track->fRow=(i&0x00ff0000)>>16;
109 track->fSector = sec;
110 // Int_t index = i&0xFFFF;
111 if (sec>=fParam->GetNInnerSector()) track->fRow += fParam->GetNRowLow();
112 track->SetClusterIndex2(track->fRow, i);
113 //track->fFirstPoint = row;
114 //if ( track->fLastPoint<row) track->fLastPoint =row;
115 // if (track->fRow<0 || track->fRow>160) {
116 // printf("problem\n");
118 if (track->fFirstPoint>track->fRow)
119 track->fFirstPoint = track->fRow;
120 if (track->fLastPoint<track->fRow)
121 track->fLastPoint = track->fRow;
124 track->fClusterPointer[track->fRow] = c;
127 Float_t angle2 = track->GetSnp()*track->GetSnp();
128 angle2 = TMath::Sqrt(angle2/(1-angle2));
130 //SET NEW Track Point
134 AliTPCTrackerPoint &point =*(track->GetTrackPoint(track->fRow));
136 point.SetSigmaY(c->GetSigmaY2()/track->fCurrentSigmaY2);
137 point.SetSigmaZ(c->GetSigmaZ2()/track->fCurrentSigmaZ2);
138 point.SetErrY(sqrt(track->fErrorY2));
139 point.SetErrZ(sqrt(track->fErrorZ2));
141 point.SetX(track->GetX());
142 point.SetY(track->GetY());
143 point.SetZ(track->GetZ());
144 point.SetAngleY(angle2);
145 point.SetAngleZ(track->GetTgl());
146 if (point.fIsShared){
147 track->fErrorY2 *= 4;
148 track->fErrorZ2 *= 4;
152 Double_t chi2 = track->GetPredictedChi2(track->fCurrentCluster);
154 track->fErrorY2 *= 1.3;
155 track->fErrorY2 += 0.01;
156 track->fErrorZ2 *= 1.3;
157 track->fErrorZ2 += 0.005;
159 if (accept>0) return 0;
160 if (track->GetNumberOfClusters()%20==0){
161 // if (track->fHelixIn){
162 // TClonesArray & larr = *(track->fHelixIn);
163 // Int_t ihelix = larr.GetEntriesFast();
164 // new(larr[ihelix]) AliHelix(*track) ;
167 track->fNoCluster =0;
168 return track->Update(c,chi2,i);
173 Int_t AliTPCtrackerMI::AcceptCluster(AliTPCseed * seed, AliTPCclusterMI * cluster, Float_t factor,
174 Float_t cory, Float_t corz)
177 // decide according desired precision to accept given
178 // cluster for tracking
179 Double_t sy2=ErrY2(seed,cluster)*cory;
180 Double_t sz2=ErrZ2(seed,cluster)*corz;
181 //sy2=ErrY2(seed,cluster)*cory;
182 //sz2=ErrZ2(seed,cluster)*cory;
184 Double_t sdistancey2 = sy2+seed->GetSigmaY2();
185 Double_t sdistancez2 = sz2+seed->GetSigmaZ2();
187 Double_t rdistancey2 = (seed->fCurrentCluster->GetY()-seed->GetY())*
188 (seed->fCurrentCluster->GetY()-seed->GetY())/sdistancey2;
189 Double_t rdistancez2 = (seed->fCurrentCluster->GetZ()-seed->GetZ())*
190 (seed->fCurrentCluster->GetZ()-seed->GetZ())/sdistancez2;
192 Double_t rdistance2 = rdistancey2+rdistancez2;
195 if (rdistance2>16) return 3;
198 if ((rdistancey2>9.*factor || rdistancez2>9.*factor) && cluster->GetType()==0)
199 return 2; //suspisiouce - will be changed
201 if ((rdistancey2>6.25*factor || rdistancez2>6.25*factor) && cluster->GetType()>0)
202 // strict cut on overlaped cluster
203 return 2; //suspisiouce - will be changed
205 if ( (rdistancey2>1.*factor || rdistancez2>6.25*factor )
206 && cluster->GetType()<0){
216 //_____________________________________________________________________________
217 AliTPCtrackerMI::AliTPCtrackerMI(const AliTPCParam *par):
218 AliTracker(), fkNIS(par->GetNInnerSector()/2), fkNOS(par->GetNOuterSector()/2)
220 //---------------------------------------------------------------------
221 // The main TPC tracker constructor
222 //---------------------------------------------------------------------
223 fInnerSec=new AliTPCSector[fkNIS];
224 fOuterSec=new AliTPCSector[fkNOS];
227 for (i=0; i<fkNIS; i++) fInnerSec[i].Setup(par,0);
228 for (i=0; i<fkNOS; i++) fOuterSec[i].Setup(par,1);
235 Int_t nrowlow = par->GetNRowLow();
236 Int_t nrowup = par->GetNRowUp();
239 for (Int_t i=0;i<nrowlow;i++){
240 fXRow[i] = par->GetPadRowRadiiLow(i);
241 fPadLength[i]= par->GetPadPitchLength(0,i);
242 fYMax[i] = fXRow[i]*TMath::Tan(0.5*par->GetInnerAngle());
246 for (Int_t i=0;i<nrowup;i++){
247 fXRow[i+nrowlow] = par->GetPadRowRadiiUp(i);
248 fPadLength[i+nrowlow] = par->GetPadPitchLength(60,i);
249 fYMax[i+nrowlow] = fXRow[i+nrowlow]*TMath::Tan(0.5*par->GetOuterAngle());
262 //_____________________________________________________________________________
263 AliTPCtrackerMI::~AliTPCtrackerMI() {
264 //------------------------------------------------------------------
265 // TPC tracker destructor
266 //------------------------------------------------------------------
275 void AliTPCtrackerMI::SetIO()
279 fInput = AliRunLoader::GetTreeR("TPC", kFALSE,AliConfig::fgkDefaultEventFolderName);
281 fOutput = AliRunLoader::GetTreeT("TPC", kTRUE,AliConfig::fgkDefaultEventFolderName);
283 AliTPCtrack *iotrack= new AliTPCtrack;
284 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
290 void AliTPCtrackerMI::SetIO(TTree * input, TTree * output, AliESD * event)
306 AliTPCtrack *iotrack= new AliTPCtrack;
307 // iotrack->fHelixIn = new TClonesArray("AliHelix");
308 //iotrack->fHelixOut = new TClonesArray("AliHelix");
309 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
312 if (output && (fDebug&2)){
313 //write the full seed information if specified in debug mode
315 fSeedTree = new TTree("Seeds","Seeds");
316 AliTPCseed * vseed = new AliTPCseed;
318 TClonesArray * arrtr = new TClonesArray("AliTPCTrackPoint",160);
319 arrtr->ExpandCreateFast(160);
320 TClonesArray * arre = new TClonesArray("AliTPCExactPoint",160);
322 vseed->fPoints = arrtr;
323 vseed->fEPoints = arre;
324 // vseed->fClusterPoints = arrcl;
325 fSeedTree->Branch("seeds","AliTPCseed",&vseed,32000,99);
328 fTreeDebug = new TTree("trackDebug","trackDebug");
329 TClonesArray * arrd = new TClonesArray("AliTPCTrackPoint2",0);
330 fTreeDebug->Branch("debug",&arrd,32000,99);
338 void AliTPCtrackerMI::FillESD(TObjArray* arr)
342 //fill esds using updated tracks
344 // write tracks to the event
345 // store index of the track
346 Int_t nseed=arr->GetEntriesFast();
347 for (Int_t i=0; i<nseed; i++) {
348 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
350 pt->PropagateTo(fParam->GetInnerRadiusLow());
351 if ( (pt->GetNumberOfClusters()>70)&& (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.55) {
353 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
354 //iotrack.SetTPCindex(i);
355 fEvent->AddTrack(&iotrack);
361 void AliTPCtrackerMI::WriteTracks(TTree * tree)
364 // write tracks from seed array to selected tree
368 AliTPCtrack *iotrack= new AliTPCtrack;
369 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
374 void AliTPCtrackerMI::WriteTracks()
377 // write tracks to the given output tree -
378 // output specified with SetIO routine
385 AliTPCtrack *iotrack= 0;
386 Int_t nseed=fSeeds->GetEntriesFast();
387 //for (Int_t i=0; i<nseed; i++) {
388 // iotrack= (AliTPCtrack*)fSeeds->UncheckedAt(i);
389 // if (iotrack) break;
391 //TBranch * br = fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
392 TBranch * br = fOutput->GetBranch("tracks");
393 br->SetAddress(&iotrack);
395 for (Int_t i=0; i<nseed; i++) {
396 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i);
398 AliTPCtrack * track = new AliTPCtrack(*pt);
401 // br->SetAddress(&iotrack);
406 //fOutput->GetDirectory()->cd();
412 //write the full seed information if specified in debug mode
414 AliTPCseed * vseed = new AliTPCseed;
416 TClonesArray * arrtr = new TClonesArray("AliTPCTrackPoint",160);
417 arrtr->ExpandCreateFast(160);
418 //TClonesArray * arrcl = new TClonesArray("AliTPCclusterMI",160);
419 //arrcl->ExpandCreateFast(160);
420 TClonesArray * arre = new TClonesArray("AliTPCExactPoint",160);
422 vseed->fPoints = arrtr;
423 vseed->fEPoints = arre;
424 // vseed->fClusterPoints = arrcl;
425 //TBranch * brseed = seedtree->Branch("seeds","AliTPCseed",&vseed,32000,99);
426 TBranch * brseed = fSeedTree->GetBranch("seeds");
428 Int_t nseed=fSeeds->GetEntriesFast();
430 for (Int_t i=0; i<nseed; i++) {
431 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i);
434 // pt->fClusterPoints = arrcl;
438 brseed->SetAddress(&vseed);
442 // pt->fClusterPoints = 0;
445 if (fTreeDebug) fTreeDebug->Write();
453 Double_t AliTPCtrackerMI::ErrY2(AliTPCseed* seed, AliTPCclusterMI * cl){
456 //seed->SetErrorY2(0.1);
458 //calculate look-up table at the beginning
459 static Bool_t ginit = kFALSE;
460 static Float_t gnoise1,gnoise2,gnoise3;
461 static Float_t ggg1[10000];
462 static Float_t ggg2[10000];
463 static Float_t ggg3[10000];
464 static Float_t glandau1[10000];
465 static Float_t glandau2[10000];
466 static Float_t glandau3[10000];
468 static Float_t gcor01[500];
469 static Float_t gcor02[500];
470 static Float_t gcorp[500];
475 for (Int_t i=1;i<500;i++){
476 Float_t rsigma = float(i)/100.;
477 gcor02[i] = TMath::Max(0.78 +TMath::Exp(7.4*(rsigma-1.2)),0.6);
478 gcor01[i] = TMath::Max(0.72 +TMath::Exp(3.36*(rsigma-1.2)),0.6);
479 gcorp[i] = TMath::Max(TMath::Power((rsigma+0.5),1.5),1.2);
483 for (Int_t i=3;i<10000;i++){
487 Float_t amp = float(i);
488 Float_t padlength =0.75;
489 gnoise1 = 0.0004/padlength;
490 Float_t nel = 0.268*amp;
491 Float_t nprim = 0.155*amp;
492 ggg1[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.001*nel/(padlength*padlength))/nel;
493 glandau1[i] = (2.+0.12*nprim)*0.5* (2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
494 if (glandau1[i]>1) glandau1[i]=1;
495 glandau1[i]*=padlength*padlength/12.;
499 gnoise2 = 0.0004/padlength;
502 ggg2[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
503 glandau2[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
504 if (glandau2[i]>1) glandau2[i]=1;
505 glandau2[i]*=padlength*padlength/12.;
510 gnoise3 = 0.0004/padlength;
513 ggg3[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
514 glandau3[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
515 if (glandau3[i]>1) glandau3[i]=1;
516 glandau3[i]*=padlength*padlength/12.;
524 Int_t amp = int(TMath::Abs(cl->GetQ()));
526 seed->SetErrorY2(1.);
530 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
531 Int_t ctype = cl->GetType();
532 Float_t padlength= GetPadPitchLength(seed->fRow);
533 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
534 angle2 = angle2/(1-angle2);
537 Int_t rsigmay = int(100.*cl->GetSigmaY2()/(seed->fCurrentSigmaY2));
540 if (fSectors==fInnerSec){
542 res = ggg1[amp]*z+glandau1[amp]*angle2;
543 if (ctype==0) res *= gcor01[rsigmay];
546 res*= gcorp[rsigmay];
552 res = ggg2[amp]*z+glandau2[amp]*angle2;
553 if (ctype==0) res *= gcor02[rsigmay];
556 res*= gcorp[rsigmay];
561 res = ggg3[amp]*z+glandau3[amp]*angle2;
562 if (ctype==0) res *= gcor02[rsigmay];
565 res*= gcorp[rsigmay];
572 res*=2.4; // overestimate error 2 times
579 seed->SetErrorY2(res);
587 Double_t AliTPCtrackerMI::ErrZ2(AliTPCseed* seed, AliTPCclusterMI * cl){
590 //seed->SetErrorY2(0.1);
592 //calculate look-up table at the beginning
593 static Bool_t ginit = kFALSE;
594 static Float_t gnoise1,gnoise2,gnoise3;
595 static Float_t ggg1[10000];
596 static Float_t ggg2[10000];
597 static Float_t ggg3[10000];
598 static Float_t glandau1[10000];
599 static Float_t glandau2[10000];
600 static Float_t glandau3[10000];
602 static Float_t gcor01[1000];
603 static Float_t gcor02[1000];
604 static Float_t gcorp[1000];
609 for (Int_t i=1;i<1000;i++){
610 Float_t rsigma = float(i)/100.;
611 gcor02[i] = TMath::Max(0.81 +TMath::Exp(6.8*(rsigma-1.2)),0.6);
612 gcor01[i] = TMath::Max(0.72 +TMath::Exp(2.04*(rsigma-1.2)),0.6);
613 gcorp[i] = TMath::Max(TMath::Power((rsigma+0.5),1.5),1.2);
617 for (Int_t i=3;i<10000;i++){
621 Float_t amp = float(i);
622 Float_t padlength =0.75;
623 gnoise1 = 0.0004/padlength;
624 Float_t nel = 0.268*amp;
625 Float_t nprim = 0.155*amp;
626 ggg1[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.001*nel/(padlength*padlength))/nel;
627 glandau1[i] = (2.+0.12*nprim)*0.5* (2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
628 if (glandau1[i]>1) glandau1[i]=1;
629 glandau1[i]*=padlength*padlength/12.;
633 gnoise2 = 0.0004/padlength;
636 ggg2[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
637 glandau2[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
638 if (glandau2[i]>1) glandau2[i]=1;
639 glandau2[i]*=padlength*padlength/12.;
644 gnoise3 = 0.0004/padlength;
647 ggg3[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
648 glandau3[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
649 if (glandau3[i]>1) glandau3[i]=1;
650 glandau3[i]*=padlength*padlength/12.;
658 Int_t amp = int(TMath::Abs(cl->GetQ()));
660 seed->SetErrorY2(1.);
664 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
665 Int_t ctype = cl->GetType();
666 Float_t padlength= GetPadPitchLength(seed->fRow);
668 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
669 // if (angle2<0.6) angle2 = 0.6;
670 angle2 = seed->GetTgl()*seed->GetTgl()*(1+angle2/(1-angle2));
673 Int_t rsigmaz = int(100.*cl->GetSigmaZ2()/(seed->fCurrentSigmaZ2));
676 if (fSectors==fInnerSec){
678 res = ggg1[amp]*z+glandau1[amp]*angle2;
679 if (ctype==0) res *= gcor01[rsigmaz];
682 res*= gcorp[rsigmaz];
688 res = ggg2[amp]*z+glandau2[amp]*angle2;
689 if (ctype==0) res *= gcor02[rsigmaz];
692 res*= gcorp[rsigmaz];
697 res = ggg3[amp]*z+glandau3[amp]*angle2;
698 if (ctype==0) res *= gcor02[rsigmaz];
701 res*= gcorp[rsigmaz];
710 if ((ctype<0) &&<70){
718 seed->SetErrorZ2(res);
725 Double_t AliTPCtrackerMI::ErrZ2(AliTPCseed* seed, AliTPCclusterMI * cl){
728 //seed->SetErrorZ2(0.1);
732 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
734 Float_t rsigmaz = cl->GetSigmaZ2()/(seed->fCurrentSigmaZ2);
735 Int_t ctype = cl->GetType();
736 Float_t amp = TMath::Abs(cl->GetQ());
741 Float_t landau=2 ; //landau fluctuation part
742 Float_t gg=2; // gg fluctuation part
743 Float_t padlength= GetPadPitchLength(seed->GetX());
745 if (fSectors==fInnerSec){
746 snoise2 = 0.0004/padlength;
749 gg = (2+0.001*nel/(padlength*padlength))/nel;
750 landau = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
751 if (landau>1) landau=1;
754 snoise2 = 0.0004/padlength;
757 gg = (2+0.0008*nel/(padlength*padlength))/nel;
758 landau = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
759 if (landau>1) landau=1;
761 Float_t sdiff = gg*fParam->GetDiffT()*fParam->GetDiffT()*z;
764 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
765 angle2 = TMath::Sqrt((1-angle2));
766 if (angle2<0.6) angle2 = 0.6;
769 Float_t angle = seed->GetTgl()/angle2;
770 Float_t angular = landau*angle*angle*padlength*padlength/12.;
771 Float_t res = sdiff + angular;
774 if ((ctype==0) && (fSectors ==fOuterSec))
775 res *= 0.81 +TMath::Exp(6.8*(rsigmaz-1.2));
777 if ((ctype==0) && (fSectors ==fInnerSec))
778 res *= 0.72 +TMath::Exp(2.04*(rsigmaz-1.2));
782 res*= TMath::Power(rsigmaz+0.5,1.5); //0.31+0.147*ctype;
788 if ((ctype<0) &&<70){
796 seed->SetErrorZ2(res);
803 void AliTPCseed::Reset(Bool_t all)
807 SetNumberOfClusters(0);
813 for (Int_t i=0;i<8;i++){
814 delete [] fTrackPoints[i];
822 for (Int_t i=0;i<200;i++) SetClusterIndex2(i,-3);
823 for (Int_t i=0;i<160;i++) fClusterPointer[i]=0;
829 void AliTPCseed::Modify(Double_t factor)
832 //------------------------------------------------------------------
833 //This function makes a track forget its history :)
834 //------------------------------------------------------------------
840 fC10*=0; fC11*=factor;
841 fC20*=0; fC21*=0; fC22*=factor;
842 fC30*=0; fC31*=0; fC32*=0; fC33*=factor;
843 fC40*=0; fC41*=0; fC42*=0; fC43*=0; fC44*=factor;
844 SetNumberOfClusters(0);
848 fCurrentSigmaY2 = 0.000005;
849 fCurrentSigmaZ2 = 0.000005;
858 Int_t AliTPCseed::GetProlongation(Double_t xk, Double_t &y, Double_t & z) const
860 //-----------------------------------------------------------------
861 // This function find proloncation of a track to a reference plane x=xk.
862 // doesn't change internal state of the track
863 //-----------------------------------------------------------------
865 Double_t x1=fX, x2=x1+(xk-x1), dx=x2-x1;
867 if (TMath::Abs(fP4*xk - fP2) >= 0.999) {
871 // Double_t y1=fP0, z1=fP1;
872 Double_t c1=fP4*x1 - fP2, r1=sqrt(1.- c1*c1);
873 Double_t c2=fP4*x2 - fP2, r2=sqrt(1.- c2*c2);
877 //y += dx*(c1+c2)/(r1+r2);
878 //z += dx*(c1+c2)/(c1*r2 + c2*r1)*fP3;
880 Double_t dy = dx*(c1+c2)/(r1+r2);
883 Double_t delta = fP4*dx*(c1+c2)/(c1*r2 + c2*r1);
885 if (TMath::Abs(delta)>0.0001){
886 dz = fP3*TMath::ASin(delta)/fP4;
888 dz = dx*fP3*(c1+c2)/(c1*r2 + c2*r1);
891 dz = fP3*AliTPCFastMath::FastAsin(delta)/fP4;
901 //_____________________________________________________________________________
902 Double_t AliTPCseed::GetPredictedChi2(const AliTPCclusterMI *c) const
904 //-----------------------------------------------------------------
905 // This function calculates a predicted chi2 increment.
906 //-----------------------------------------------------------------
907 //Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
908 Double_t r00=fErrorY2, r01=0., r11=fErrorZ2;
909 r00+=fC00; r01+=fC10; r11+=fC11;
911 Double_t det=r00*r11 - r01*r01;
912 if (TMath::Abs(det) < 1.e-10) {
913 Int_t n=GetNumberOfClusters();
914 if (n>4) cerr<<n<<" AliKalmanTrack warning: Singular matrix !\n";
917 Double_t tmp=r00; r00=r11; r11=tmp; r01=-r01;
919 Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
921 return (dy*r00*dy + 2*r01*dy*dz + dz*r11*dz)/det;
925 //_________________________________________________________________________________________
928 Int_t AliTPCseed::Compare(const TObject *o) const {
929 //-----------------------------------------------------------------
930 // This function compares tracks according to the sector - for given sector according z
931 //-----------------------------------------------------------------
932 AliTPCseed *t=(AliTPCseed*)o;
935 if (t->fRelativeSector>fRelativeSector) return -1;
936 if (t->fRelativeSector<fRelativeSector) return 1;
937 Double_t z2 = t->GetZ();
938 Double_t z1 = GetZ();
940 if (z2<z1) return -1;
945 f2 = 1-20*TMath::Sqrt(t->fC44)/(TMath::Abs(t->GetC())+0.0066);
946 if (t->fBConstrain) f2=1.2;
949 f1 = 1-20*TMath::Sqrt(fC44)/(TMath::Abs(GetC())+0.0066);
951 if (fBConstrain) f1=1.2;
953 if (t->GetNumberOfClusters()*f2 <GetNumberOfClusters()*f1) return -1;
958 void AliTPCtrackerMI::RotateToLocal(AliTPCseed *seed)
960 //rotate to track "local coordinata
961 Float_t x = seed->GetX();
962 Float_t y = seed->GetY();
963 Float_t ymax = x*TMath::Tan(0.5*fSectors->GetAlpha());
966 seed->fRelativeSector= (seed->fRelativeSector+1) % fN;
967 if (!seed->Rotate(fSectors->GetAlpha()))
969 } else if (y <-ymax) {
970 seed->fRelativeSector= (seed->fRelativeSector-1+fN) % fN;
971 if (!seed->Rotate(-fSectors->GetAlpha()))
980 //_____________________________________________________________________________
981 Int_t AliTPCseed::Update(const AliTPCclusterMI *c, Double_t chisq, UInt_t /*index*/) {
982 //-----------------------------------------------------------------
983 // This function associates a cluster with this track.
984 //-----------------------------------------------------------------
985 Double_t r00=fErrorY2, r01=0., r11=fErrorZ2;
987 r00+=fC00; r01+=fC10; r11+=fC11;
988 Double_t det=r00*r11 - r01*r01;
989 Double_t tmp=r00; r00=r11/det; r11=tmp/det; r01=-r01/det;
991 Double_t k00=fC00*r00+fC10*r01, k01=fC00*r01+fC10*r11;
992 Double_t k10=fC10*r00+fC11*r01, k11=fC10*r01+fC11*r11;
993 Double_t k20=fC20*r00+fC21*r01, k21=fC20*r01+fC21*r11;
994 Double_t k30=fC30*r00+fC31*r01, k31=fC30*r01+fC31*r11;
995 Double_t k40=fC40*r00+fC41*r01, k41=fC40*r01+fC41*r11;
997 Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
998 Double_t cur=fP4 + k40*dy + k41*dz, eta=fP2 + k20*dy + k21*dz;
999 if (TMath::Abs(cur*fX-eta) >= 0.9) {
1003 fP0 += k00*dy + k01*dz;
1004 fP1 += k10*dy + k11*dz;
1006 fP3 += k30*dy + k31*dz;
1009 Double_t c01=fC10, c02=fC20, c03=fC30, c04=fC40;
1010 Double_t c12=fC21, c13=fC31, c14=fC41;
1012 fC00-=k00*fC00+k01*fC10; fC10-=k00*c01+k01*fC11;
1013 fC20-=k00*c02+k01*c12; fC30-=k00*c03+k01*c13;
1014 fC40-=k00*c04+k01*c14;
1016 fC11-=k10*c01+k11*fC11;
1017 fC21-=k10*c02+k11*c12; fC31-=k10*c03+k11*c13;
1018 fC41-=k10*c04+k11*c14;
1020 fC22-=k20*c02+k21*c12; fC32-=k20*c03+k21*c13;
1021 fC42-=k20*c04+k21*c14;
1023 fC33-=k30*c03+k31*c13;
1024 fC43-=k40*c03+k41*c13;
1026 fC44-=k40*c04+k41*c14;
1028 Int_t n=GetNumberOfClusters();
1030 SetNumberOfClusters(n+1);
1031 SetChi2(GetChi2()+chisq);
1038 //_____________________________________________________________________________
1039 Double_t AliTPCtrackerMI::F1old(Double_t x1,Double_t y1,
1040 Double_t x2,Double_t y2,
1041 Double_t x3,Double_t y3)
1043 //-----------------------------------------------------------------
1044 // Initial approximation of the track curvature
1045 //-----------------------------------------------------------------
1046 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
1047 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
1048 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
1049 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
1050 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
1052 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
1053 if ( xr*xr+yr*yr<=0.00000000000001) return 100;
1054 return -xr*yr/sqrt(xr*xr+yr*yr);
1059 //_____________________________________________________________________________
1060 Double_t AliTPCtrackerMI::F1(Double_t x1,Double_t y1,
1061 Double_t x2,Double_t y2,
1062 Double_t x3,Double_t y3)
1064 //-----------------------------------------------------------------
1065 // Initial approximation of the track curvature
1066 //-----------------------------------------------------------------
1072 Double_t det = x3*y2-x2*y3;
1077 Double_t u = 0.5* (x2*(x2-x3)+y2*(y2-y3))/det;
1078 Double_t x0 = x3*0.5-y3*u;
1079 Double_t y0 = y3*0.5+x3*u;
1080 Double_t c2 = 1/TMath::Sqrt(x0*x0+y0*y0);
1086 Double_t AliTPCtrackerMI::F2(Double_t x1,Double_t y1,
1087 Double_t x2,Double_t y2,
1088 Double_t x3,Double_t y3)
1090 //-----------------------------------------------------------------
1091 // Initial approximation of the track curvature
1092 //-----------------------------------------------------------------
1098 Double_t det = x3*y2-x2*y3;
1103 Double_t u = 0.5* (x2*(x2-x3)+y2*(y2-y3))/det;
1104 Double_t x0 = x3*0.5-y3*u;
1105 Double_t y0 = y3*0.5+x3*u;
1106 Double_t c2 = 1/TMath::Sqrt(x0*x0+y0*y0);
1115 //_____________________________________________________________________________
1116 Double_t AliTPCtrackerMI::F2old(Double_t x1,Double_t y1,
1117 Double_t x2,Double_t y2,
1118 Double_t x3,Double_t y3)
1120 //-----------------------------------------------------------------
1121 // Initial approximation of the track curvature times center of curvature
1122 //-----------------------------------------------------------------
1123 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
1124 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
1125 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
1126 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
1127 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
1129 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
1131 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
1134 //_____________________________________________________________________________
1135 Double_t AliTPCtrackerMI::F3(Double_t x1,Double_t y1,
1136 Double_t x2,Double_t y2,
1137 Double_t z1,Double_t z2)
1139 //-----------------------------------------------------------------
1140 // Initial approximation of the tangent of the track dip angle
1141 //-----------------------------------------------------------------
1142 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
1146 Double_t AliTPCtrackerMI::F3n(Double_t x1,Double_t y1,
1147 Double_t x2,Double_t y2,
1148 Double_t z1,Double_t z2, Double_t c)
1150 //-----------------------------------------------------------------
1151 // Initial approximation of the tangent of the track dip angle
1152 //-----------------------------------------------------------------
1156 //angle1 = (z1-z2)*c/(TMath::ASin(c*x1-ni)-TMath::ASin(c*x2-ni));
1158 Double_t d = TMath::Sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
1159 if (TMath::Abs(d*c*0.5)>1) return 0;
1160 // Double_t angle2 = TMath::ASin(d*c*0.5);
1161 // Double_t angle2 = AliTPCFastMath::FastAsin(d*c*0.5);
1162 Double_t angle2 = (d*c*0.5>0.1)? TMath::ASin(d*c*0.5): AliTPCFastMath::FastAsin(d*c*0.5);
1164 angle2 = (z1-z2)*c/(angle2*2.);
1168 Bool_t AliTPCtrackerMI::GetProlongation(Double_t x1, Double_t x2, Double_t x[5], Double_t &y, Double_t &z)
1169 {//-----------------------------------------------------------------
1170 // This function find proloncation of a track to a reference plane x=x2.
1171 //-----------------------------------------------------------------
1175 if (TMath::Abs(x[4]*x1 - x[2]) >= 0.999) {
1179 Double_t c1=x[4]*x1 - x[2], r1=sqrt(1.- c1*c1);
1180 Double_t c2=x[4]*x2 - x[2], r2=sqrt(1.- c2*c2);
1184 Double_t dy = dx*(c1+c2)/(r1+r2);
1187 Double_t delta = x[4]*dx*(c1+c2)/(c1*r2 + c2*r1);
1189 if (TMath::Abs(delta)>0.01){
1190 dz = x[3]*TMath::ASin(delta)/x[4];
1192 dz = x[3]*AliTPCFastMath::FastAsin(delta)/x[4];
1195 //dz = x[3]*AliTPCFastMath::FastAsin(delta)/x[4];
1203 Int_t AliTPCtrackerMI::LoadClusters (TTree *tree)
1208 return LoadClusters();
1211 Int_t AliTPCtrackerMI::LoadClusters()
1214 // load clusters to the memory
1215 AliTPCClustersRow *clrow= new AliTPCClustersRow;
1216 clrow->SetClass("AliTPCclusterMI");
1218 clrow->GetArray()->ExpandCreateFast(10000);
1220 // TTree * tree = fClustersArray.GetTree();
1222 TTree * tree = fInput;
1223 TBranch * br = tree->GetBranch("Segment");
1224 br->SetAddress(&clrow);
1226 Int_t j=Int_t(tree->GetEntries());
1227 for (Int_t i=0; i<j; i++) {
1231 fParam->AdjustSectorRow(clrow->GetID(),sec,row);
1233 AliTPCRow * tpcrow=0;
1236 tpcrow = &(fInnerSec[sec%fkNIS][row]);
1240 tpcrow = &(fOuterSec[(sec-fkNIS*2)%fkNOS][row]);
1241 left = (sec-fkNIS*2)/fkNOS;
1244 tpcrow->fN1 = clrow->GetArray()->GetEntriesFast();
1245 tpcrow->fClusters1 = new AliTPCclusterMI[tpcrow->fN1];
1246 for (Int_t i=0;i<tpcrow->fN1;i++)
1247 tpcrow->fClusters1[i] = *(AliTPCclusterMI*)(clrow->GetArray()->At(i));
1250 tpcrow->fN2 = clrow->GetArray()->GetEntriesFast();
1251 tpcrow->fClusters2 = new AliTPCclusterMI[tpcrow->fN2];
1252 for (Int_t i=0;i<tpcrow->fN2;i++)
1253 tpcrow->fClusters2[i] = *(AliTPCclusterMI*)(clrow->GetArray()->At(i));
1264 void AliTPCtrackerMI::UnloadClusters()
1267 // unload clusters from the memory
1269 Int_t nrows = fOuterSec->GetNRows();
1270 for (Int_t sec = 0;sec<fkNOS;sec++)
1271 for (Int_t row = 0;row<nrows;row++){
1272 AliTPCRow* tpcrow = &(fOuterSec[sec%fkNOS][row]);
1274 // if (tpcrow->fClusters1) delete []tpcrow->fClusters1;
1275 // if (tpcrow->fClusters2) delete []tpcrow->fClusters2;
1277 tpcrow->ResetClusters();
1280 nrows = fInnerSec->GetNRows();
1281 for (Int_t sec = 0;sec<fkNIS;sec++)
1282 for (Int_t row = 0;row<nrows;row++){
1283 AliTPCRow* tpcrow = &(fInnerSec[sec%fkNIS][row]);
1285 // if (tpcrow->fClusters1) delete []tpcrow->fClusters1;
1286 //if (tpcrow->fClusters2) delete []tpcrow->fClusters2;
1288 tpcrow->ResetClusters();
1295 //_____________________________________________________________________________
1296 Int_t AliTPCtrackerMI::LoadOuterSectors() {
1297 //-----------------------------------------------------------------
1298 // This function fills outer TPC sectors with clusters.
1299 //-----------------------------------------------------------------
1300 Int_t nrows = fOuterSec->GetNRows();
1302 for (Int_t sec = 0;sec<fkNOS;sec++)
1303 for (Int_t row = 0;row<nrows;row++){
1304 AliTPCRow* tpcrow = &(fOuterSec[sec%fkNOS][row]);
1305 Int_t sec2 = sec+2*fkNIS;
1307 Int_t ncl = tpcrow->fN1;
1309 AliTPCclusterMI *c= &(tpcrow->fClusters1[ncl]);
1310 index=(((sec2<<8)+row)<<16)+ncl;
1311 tpcrow->InsertCluster(c,index);
1316 AliTPCclusterMI *c= &(tpcrow->fClusters2[ncl]);
1317 index=((((sec2+fkNOS)<<8)+row)<<16)+ncl;
1318 tpcrow->InsertCluster(c,index);
1321 // write indexes for fast acces
1323 for (Int_t i=0;i<510;i++)
1324 tpcrow->fFastCluster[i]=-1;
1325 for (Int_t i=0;i<tpcrow->GetN();i++){
1326 Int_t zi = Int_t((*tpcrow)[i]->GetZ()+255.);
1327 tpcrow->fFastCluster[zi]=i; // write index
1330 for (Int_t i=0;i<510;i++){
1331 if (tpcrow->fFastCluster[i]<0)
1332 tpcrow->fFastCluster[i] = last;
1334 last = tpcrow->fFastCluster[i];
1343 //_____________________________________________________________________________
1344 Int_t AliTPCtrackerMI::LoadInnerSectors() {
1345 //-----------------------------------------------------------------
1346 // This function fills inner TPC sectors with clusters.
1347 //-----------------------------------------------------------------
1348 Int_t nrows = fInnerSec->GetNRows();
1350 for (Int_t sec = 0;sec<fkNIS;sec++)
1351 for (Int_t row = 0;row<nrows;row++){
1352 AliTPCRow* tpcrow = &(fInnerSec[sec%fkNIS][row]);
1355 Int_t ncl = tpcrow->fN1;
1357 AliTPCclusterMI *c= &(tpcrow->fClusters1[ncl]);
1358 index=(((sec<<8)+row)<<16)+ncl;
1359 tpcrow->InsertCluster(c,index);
1364 AliTPCclusterMI *c= &(tpcrow->fClusters2[ncl]);
1365 index=((((sec+fkNIS)<<8)+row)<<16)+ncl;
1366 tpcrow->InsertCluster(c,index);
1369 // write indexes for fast acces
1371 for (Int_t i=0;i<510;i++)
1372 tpcrow->fFastCluster[i]=-1;
1373 for (Int_t i=0;i<tpcrow->GetN();i++){
1374 Int_t zi = Int_t((*tpcrow)[i]->GetZ()+255.);
1375 tpcrow->fFastCluster[zi]=i; // write index
1378 for (Int_t i=0;i<510;i++){
1379 if (tpcrow->fFastCluster[i]<0)
1380 tpcrow->fFastCluster[i] = last;
1382 last = tpcrow->fFastCluster[i];
1394 //_________________________________________________________________________
1395 AliTPCclusterMI *AliTPCtrackerMI::GetClusterMI(Int_t index) const {
1396 //--------------------------------------------------------------------
1397 // Return pointer to a given cluster
1398 //--------------------------------------------------------------------
1399 Int_t sec=(index&0xff000000)>>24;
1400 Int_t row=(index&0x00ff0000)>>16;
1401 Int_t ncl=(index&0x00007fff)>>00;
1403 const AliTPCRow * tpcrow=0;
1404 AliTPCclusterMI * clrow =0;
1406 tpcrow = &(fInnerSec[sec%fkNIS][row]);
1408 clrow = tpcrow->fClusters1;
1410 clrow = tpcrow->fClusters2;
1413 tpcrow = &(fOuterSec[(sec-fkNIS*2)%fkNOS][row]);
1414 if (sec-2*fkNIS<fkNOS)
1415 clrow = tpcrow->fClusters1;
1417 clrow = tpcrow->fClusters2;
1419 if (tpcrow==0) return 0;
1420 if (tpcrow->GetN()<=ncl) return 0;
1421 // return (AliTPCclusterMI*)(*tpcrow)[ncl];
1422 return &(clrow[ncl]);
1428 Int_t AliTPCtrackerMI::FollowToNext(AliTPCseed& t, Int_t nr) {
1429 //-----------------------------------------------------------------
1430 // This function tries to find a track prolongation to next pad row
1431 //-----------------------------------------------------------------
1433 Double_t x= GetXrow(nr), ymax=GetMaxY(nr);
1434 AliTPCclusterMI *cl=0;
1435 Int_t tpcindex= t.GetClusterIndex2(nr);
1437 // update current shape info every 5 pad-row
1438 // if ( (nr%5==0) || t.GetNumberOfClusters()<2 || (t.fCurrentSigmaY2<0.0001) ){
1442 if (fIteration>0 && tpcindex>=-1){ //if we have already clusters
1444 if (tpcindex==-1) return 0; //track in dead zone
1446 cl = t.fClusterPointer[nr];
1447 if ( (cl==0) && (index>0)) cl = GetClusterMI(tpcindex);
1448 t.fCurrentClusterIndex1 = tpcindex;
1451 Int_t relativesector = ((tpcindex&0xff000000)>>24)%18; // if previously accepted cluster in different sector
1452 Float_t angle = relativesector*fSectors->GetAlpha()+fSectors->GetAlphaShift();
1454 if (angle<-TMath::Pi()) angle += 2*TMath::Pi();
1455 if (angle>=TMath::Pi()) angle -= 2*TMath::Pi();
1457 if (TMath::Abs(angle-t.GetAlpha())>0.001){
1458 Double_t rotation = angle-t.GetAlpha();
1459 t.fRelativeSector= relativesector;
1464 t.fCurrentCluster = cl;
1466 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1467 if ((tpcindex&0x8000)==0) accept =0;
1469 //if founded cluster is acceptible
1470 if (cl->IsUsed(11)) { // id cluster is shared inrease uncertainty
1477 UpdateTrack(&t,accept);
1482 if (fIteration>1) return 0; // not look for new cluster during refitting
1485 if (TMath::Abs(t.GetSnp())>0.95 || TMath::Abs(x*t.GetC()-t.GetEta())>0.95) return 0;
1486 Double_t y=t.GetYat(x);
1487 if (TMath::Abs(y)>ymax){
1489 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1490 if (!t.Rotate(fSectors->GetAlpha()))
1492 } else if (y <-ymax) {
1493 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1494 if (!t.Rotate(-fSectors->GetAlpha()))
1500 if (!t.PropagateTo(x)) {
1501 if (fIteration==0) t.fRemoval = 10;
1505 Double_t z=t.GetZ();
1507 const AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1508 if ( (t.GetSigmaY2()<0) || t.GetSigmaZ2()<0) return 0;
1510 Double_t roadz = 1.;
1512 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1514 t.SetClusterIndex2(nr,-1);
1519 if (TMath::Abs(z)<(1.05*x+10)) t.fNFoundable++;
1525 // cl = krow.FindNearest2(y+10.,z,roady,roadz,index);
1526 cl = krow.FindNearest2(y,z,roady,roadz,index);
1527 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1530 t.fCurrentCluster = cl;
1532 if (fIteration==2&&cl->IsUsed(10)) return 0;
1533 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1534 if (fIteration==2&&cl->IsUsed(11)) {
1541 if (t.fCurrentCluster->IsUsed(10)){
1546 if (t.fNShared>0.7*t.GetNumberOfClusters()) {
1552 if (accept<3) UpdateTrack(&t,accept);
1555 if ( fIteration==0 && t.fNFoundable*0.5 > t.GetNumberOfClusters()) t.fRemoval=10;
1561 Int_t AliTPCtrackerMI::FollowToNextFast(AliTPCseed& t, Int_t nr) {
1562 //-----------------------------------------------------------------
1563 // This function tries to find a track prolongation to next pad row
1564 //-----------------------------------------------------------------
1566 Double_t x= GetXrow(nr), ymax=GetMaxY(nr);
1568 if (!t.GetProlongation(x,y,z)) {
1574 if (TMath::Abs(y)>ymax){
1577 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1578 if (!t.Rotate(fSectors->GetAlpha()))
1580 } else if (y <-ymax) {
1581 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1582 if (!t.Rotate(-fSectors->GetAlpha()))
1585 if (!t.PropagateTo(x)) {
1588 t.GetProlongation(x,y,z);
1591 // update current shape info every 3 pad-row
1592 if ( (nr%6==0) || t.GetNumberOfClusters()<2 || (t.fCurrentSigmaY2<0.0001) ){
1593 // t.fCurrentSigmaY = GetSigmaY(&t);
1594 //t.fCurrentSigmaZ = GetSigmaZ(&t);
1598 AliTPCclusterMI *cl=0;
1603 const AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1604 if ( (t.GetSigmaY2()<0) || t.GetSigmaZ2()<0) return 0;
1606 Double_t roadz = 1.;
1609 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1611 t.SetClusterIndex2(row,-1);
1616 if (TMath::Abs(z)>(1.05*x+10)) t.SetClusterIndex2(row,-1);
1620 if ((cl==0)&&(krow)) {
1621 // cl = krow.FindNearest2(y+10,z,roady,roadz,index);
1622 cl = krow.FindNearest2(y,z,roady,roadz,index);
1624 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1628 t.fCurrentCluster = cl;
1629 // Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1631 t.SetClusterIndex2(row,index);
1632 t.fClusterPointer[row] = cl;
1640 Int_t AliTPCtrackerMI::UpdateClusters(AliTPCseed& t, Int_t nr) {
1641 //-----------------------------------------------------------------
1642 // This function tries to find a track prolongation to next pad row
1643 //-----------------------------------------------------------------
1644 t.fCurrentCluster = 0;
1645 t.fCurrentClusterIndex1 = 0;
1647 Double_t xt=t.GetX();
1648 Int_t row = GetRowNumber(xt)-1;
1649 Double_t ymax= GetMaxY(nr);
1651 if (row < nr) return 1; // don't prolongate if not information until now -
1652 if (TMath::Abs(t.GetSnp())>0.9 && t.GetNumberOfClusters()>40. && fIteration!=2) {
1654 return 0; // not prolongate strongly inclined tracks
1656 if (TMath::Abs(t.GetSnp())>0.95) {
1658 return 0; // not prolongate strongly inclined tracks
1661 Double_t x= GetXrow(nr);
1663 //t.PropagateTo(x+0.02);
1664 //t.PropagateTo(x+0.01);
1665 if (!t.PropagateTo(x)){
1672 if (TMath::Abs(y)>ymax){
1674 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1675 if (!t.Rotate(fSectors->GetAlpha()))
1677 } else if (y <-ymax) {
1678 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1679 if (!t.Rotate(-fSectors->GetAlpha()))
1682 // if (!t.PropagateTo(x)){
1690 AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1692 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1694 t.SetClusterIndex2(nr,-1);
1699 if (TMath::Abs(t.GetZ())<(1.05*t.GetX()+10)) t.fNFoundable++;
1705 if ( (nr%6==0) || t.GetNumberOfClusters()<2){
1706 // t.fCurrentSigmaY = GetSigmaY(&t);
1707 //t.fCurrentSigmaZ = GetSigmaZ(&t);
1711 AliTPCclusterMI *cl=0;
1714 Double_t roady = 1.;
1715 Double_t roadz = 1.;
1719 index = t.GetClusterIndex2(nr);
1720 if ( (index>0) && (index&0x8000)==0){
1721 cl = t.fClusterPointer[nr];
1722 if ( (cl==0) && (index>0)) cl = GetClusterMI(index);
1723 t.fCurrentClusterIndex1 = index;
1725 t.fCurrentCluster = cl;
1732 //cl = krow.FindNearest2(y+10,z,roady,roadz,index);
1733 cl = krow.FindNearest2(y,z,roady,roadz,index);
1736 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1737 t.fCurrentCluster = cl;
1743 Int_t AliTPCtrackerMI::FollowToNextCluster(AliTPCseed & t, Int_t nr) {
1744 //-----------------------------------------------------------------
1745 // This function tries to find a track prolongation to next pad row
1746 //-----------------------------------------------------------------
1748 //update error according neighborhoud
1750 if (t.fCurrentCluster) {
1752 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1754 if (t.fCurrentCluster->IsUsed(10)){
1760 if (t.fNShared>0.7*t.GetNumberOfClusters()) {
1765 if (fIteration>0) accept = 0;
1766 if (accept<3) UpdateTrack(&t,accept);
1770 if ( ( (t.GetSigmaY2()+t.GetSigmaZ2())>0.16)&& t.GetNumberOfClusters()>18) t.fRemoval=10;
1771 if ( t.GetChi2()/t.GetNumberOfClusters()>6 &&t.GetNumberOfClusters()>18) t.fRemoval=10;
1773 if (( (t.fNFoundable*0.5 > t.GetNumberOfClusters()) || t.fNoCluster>15)) t.fRemoval=10;
1781 //_____________________________________________________________________________
1782 Int_t AliTPCtrackerMI::FollowProlongation(AliTPCseed& t, Int_t rf, Int_t step) {
1783 //-----------------------------------------------------------------
1784 // This function tries to find a track prolongation.
1785 //-----------------------------------------------------------------
1786 Double_t xt=t.GetX();
1788 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1789 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1790 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1792 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1794 Int_t first = GetRowNumber(xt)-1;
1795 for (Int_t nr= first; nr>=rf; nr-=step) {
1796 if (nr<fInnerSec->GetNRows())
1797 fSectors = fInnerSec;
1799 fSectors = fOuterSec;
1800 if (FollowToNext(t,nr)==0)
1809 //_____________________________________________________________________________
1810 Int_t AliTPCtrackerMI::FollowProlongationFast(AliTPCseed& t, Int_t rf, Int_t step) {
1811 //-----------------------------------------------------------------
1812 // This function tries to find a track prolongation.
1813 //-----------------------------------------------------------------
1814 Double_t xt=t.GetX();
1816 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1817 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1818 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1819 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1821 for (Int_t nr=GetRowNumber(xt)-1; nr>=rf; nr-=step) {
1823 if (FollowToNextFast(t,nr)==0)
1824 if (!t.IsActive()) return 0;
1834 Int_t AliTPCtrackerMI::FollowBackProlongation(AliTPCseed& t, Int_t rf) {
1835 //-----------------------------------------------------------------
1836 // This function tries to find a track prolongation.
1837 //-----------------------------------------------------------------
1838 // Double_t xt=t.GetX();
1840 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1841 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1842 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1843 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1845 Int_t first = t.fFirstPoint;
1847 if (first<0) first=0;
1848 for (Int_t nr=first; nr<=rf; nr++) {
1849 //if ( (t.GetSnp()<0.9))
1850 if (nr<fInnerSec->GetNRows())
1851 fSectors = fInnerSec;
1853 fSectors = fOuterSec;
1863 Float_t AliTPCtrackerMI::OverlapFactor(AliTPCseed * s1, AliTPCseed * s2, Int_t &sum1, Int_t & sum2)
1871 Float_t dz2 =(s1->GetZ() - s2->GetZ());
1874 Float_t dy2 =TMath::Abs((s1->GetY() - s2->GetY()));
1876 Float_t distance = TMath::Sqrt(dz2+dy2);
1877 if (distance>4.) return 0; // if there are far away - not overlap - to reduce combinatorics
1880 Int_t firstpoint = TMath::Min(s1->fFirstPoint,s2->fFirstPoint);
1881 Int_t lastpoint = TMath::Max(s1->fLastPoint,s2->fLastPoint);
1886 if (firstpoint>lastpoint) {
1887 firstpoint =lastpoint;
1892 for (Int_t i=firstpoint-1;i<lastpoint+1;i++){
1893 if (s1->GetClusterIndex2(i)>0) sum1++;
1894 if (s2->GetClusterIndex2(i)>0) sum2++;
1895 if (s1->GetClusterIndex2(i)==s2->GetClusterIndex2(i) && s1->GetClusterIndex2(i)>0) {
1899 if (sum<5) return 0;
1901 Float_t summin = TMath::Min(sum1+1,sum2+1);
1902 Float_t ratio = (sum+1)/Float_t(summin);
1906 void AliTPCtrackerMI::SignShared(AliTPCseed * s1, AliTPCseed * s2)
1910 if (TMath::Abs(s1->GetC()-s2->GetC())>0.004) return;
1911 if (TMath::Abs(s1->GetTgl()-s2->GetTgl())>0.6) return;
1913 Float_t dz2 =(s1->GetZ() - s2->GetZ());
1915 Float_t dy2 =(s1->GetY() - s2->GetY());
1917 Float_t distance = dz2+dy2;
1918 if (distance>325.) return ; // if there are far away - not overlap - to reduce combinatorics
1923 Int_t firstpoint = TMath::Max(s1->fFirstPoint,s2->fFirstPoint);
1924 Int_t lastpoint = TMath::Min(s1->fLastPoint,s2->fLastPoint);
1926 if (firstpoint>=lastpoint-5) return;;
1928 for (Int_t i=firstpoint;i<lastpoint;i++){
1929 // if ( (s1->GetClusterIndex2(i)&0xFFFF8FFF)==(s2->GetClusterIndex2(i)&0xFFFF8FFF) && s1->GetClusterIndex2(i)>0) {
1930 if ( (s1->GetClusterIndex2(i))==(s2->GetClusterIndex2(i)) && s1->GetClusterIndex2(i)>0) {
1937 for (Int_t i=firstpoint;i<lastpoint;i++){
1938 // if ( (s1->GetClusterIndex2(i)&0xFFFF8FFF)==(s2->GetClusterIndex2(i)&0xFFFF8FFF) && s1->GetClusterIndex2(i)>0) {
1939 if ( (s1->GetClusterIndex2(i))==(s2->GetClusterIndex2(i)) && s1->GetClusterIndex2(i)>0) {
1940 AliTPCTrackerPoint *p1 = s1->GetTrackPoint(i);
1941 AliTPCTrackerPoint *p2 = s2->GetTrackPoint(i);;
1942 if (s1->IsActive()&&s2->IsActive()){
1943 p1->fIsShared = kTRUE;
1944 p2->fIsShared = kTRUE;
1951 for (Int_t i=0;i<4;i++){
1952 if (s1->fOverlapLabels[3*i]==0){
1953 s1->fOverlapLabels[3*i] = s2->GetLabel();
1954 s1->fOverlapLabels[3*i+1] = sumshared;
1955 s1->fOverlapLabels[3*i+2] = s2->GetUniqueID();
1959 for (Int_t i=0;i<4;i++){
1960 if (s2->fOverlapLabels[3*i]==0){
1961 s2->fOverlapLabels[3*i] = s1->GetLabel();
1962 s2->fOverlapLabels[3*i+1] = sumshared;
1963 s2->fOverlapLabels[3*i+2] = s1->GetUniqueID();
1971 void AliTPCtrackerMI::SignShared(TObjArray * arr)
1974 //sort trackss according sectors
1976 for (Int_t i=0; i<arr->GetEntriesFast(); i++) {
1977 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
1979 //if (pt) RotateToLocal(pt);
1983 arr->Sort(); // sorting according z
1984 arr->Expand(arr->GetEntries());
1987 Int_t nseed=arr->GetEntriesFast();
1988 for (Int_t i=0; i<nseed; i++) {
1989 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
1991 for (Int_t j=0;j<=12;j++){
1992 pt->fOverlapLabels[j] =0;
1995 for (Int_t i=0; i<nseed; i++) {
1996 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
1998 if (pt->fRemoval>10) continue;
1999 for (Int_t j=i+1; j<nseed; j++){
2000 AliTPCseed *pt2=(AliTPCseed*)arr->UncheckedAt(j);
2002 if (pt2->fRemoval<=10) {
2003 if ( TMath::Abs(pt->fRelativeSector-pt2->fRelativeSector)>0) break;
2010 void AliTPCtrackerMI::RemoveDouble(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t removalindex)
2013 //sort trackss according sectors
2016 Info("RemoveDouble","Number of tracks before double removal- %d\n",arr->GetEntries());
2019 for (Int_t i=0; i<arr->GetEntriesFast(); i++) {
2020 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2025 arr->Sort(); // sorting according z
2026 arr->Expand(arr->GetEntries());
2028 //reset overlap labels
2030 Int_t nseed=arr->GetEntriesFast();
2031 for (Int_t i=0; i<nseed; i++) {
2032 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2035 for (Int_t j=0;j<=12;j++){
2036 pt->fOverlapLabels[j] =0;
2040 //sign shared tracks
2041 for (Int_t i=0; i<nseed; i++) {
2042 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2044 if (pt->fRemoval>10) continue;
2045 Float_t deltac = pt->GetC()*0.1;
2046 for (Int_t j=i+1; j<nseed; j++){
2047 AliTPCseed *pt2=(AliTPCseed*)arr->UncheckedAt(j);
2049 if (pt2->fRemoval<=10) {
2050 if ( TMath::Abs(pt->fRelativeSector-pt2->fRelativeSector)>0) break;
2051 if (TMath::Abs(pt->GetC() -pt2->GetC())>deltac) continue;
2052 if (TMath::Abs(pt->GetTgl()-pt2->GetTgl())>0.05) continue;
2059 // remove highly shared tracks
2060 for (Int_t i=0; i<nseed; i++) {
2061 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2063 if (pt->fRemoval>10) continue;
2066 for (Int_t j=0;j<4;j++){
2067 sumshared = pt->fOverlapLabels[j*3+1];
2069 Float_t factor = factor1;
2070 if (pt->fRemoval>0) factor = factor2;
2071 if (sumshared/pt->GetNumberOfClusters()>factor){
2072 for (Int_t j=0;j<4;j++){
2073 if (pt->fOverlapLabels[3*j]==0) continue;
2074 if (pt->fOverlapLabels[3*j+1]<5) continue;
2075 if (pt->fRemoval==removalindex) continue;
2076 AliTPCseed * pt2 = (AliTPCseed*)arr->UncheckedAt(pt->fOverlapLabels[3*j+2]);
2078 if (pt2->GetSigma2C()<pt->GetSigma2C()){
2079 // pt->fRemoval = removalindex;
2080 delete arr->RemoveAt(i);
2088 Info("RemoveDouble","Number of tracks after double removal- %d\n",arr->GetEntries());
2097 void AliTPCtrackerMI::SortTracks(TObjArray * arr, Int_t mode) const
2100 //sort tracks in array according mode criteria
2101 Int_t nseed = arr->GetEntriesFast();
2102 for (Int_t i=0; i<nseed; i++) {
2103 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2113 void AliTPCtrackerMI::RemoveUsed(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t removalindex)
2116 //Loop over all tracks and remove "overlaps"
2119 Int_t nseed = arr->GetEntriesFast();
2122 for (Int_t i=0; i<nseed; i++) {
2123 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2125 delete arr->RemoveAt(i);
2129 pt->fBSigned = kFALSE;
2133 nseed = arr->GetEntriesFast();
2140 for (Int_t i=0; i<nseed; i++) {
2141 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2145 Int_t found,foundable,shared;
2147 pt->GetClusterStatistic(0,160,found, foundable,shared,kFALSE);
2149 pt->GetClusterStatistic(0,160,found, foundable,shared,kTRUE);
2151 Double_t factor = factor2;
2152 if (pt->fBConstrain) factor = factor1;
2154 if ((Float_t(shared)/Float_t(found))>factor){
2155 pt->Desactivate(removalindex);
2160 for (Int_t i=0; i<160; i++) {
2161 Int_t index=pt->GetClusterIndex2(i);
2162 if (index<0 || index&0x8000 ) continue;
2163 AliTPCclusterMI *c= pt->fClusterPointer[i];
2165 // if (!c->IsUsed(10)) c->Use(10);
2166 //if (pt->IsActive())
2175 Info("RemoveUsed","\n*****\nNumber of good tracks after shared removal\t%d\n",fNtracks);
2179 void AliTPCtrackerMI::UnsignClusters()
2182 // loop over all clusters and unsign them
2185 for (Int_t sec=0;sec<fkNIS;sec++){
2186 for (Int_t row=0;row<fInnerSec->GetNRows();row++){
2187 AliTPCclusterMI *cl = fInnerSec[sec][row].fClusters1;
2188 for (Int_t icl =0;icl< fInnerSec[sec][row].fN1;icl++)
2189 // if (cl[icl].IsUsed(10))
2191 cl = fInnerSec[sec][row].fClusters2;
2192 for (Int_t icl =0;icl< fInnerSec[sec][row].fN2;icl++)
2193 //if (cl[icl].IsUsed(10))
2198 for (Int_t sec=0;sec<fkNOS;sec++){
2199 for (Int_t row=0;row<fOuterSec->GetNRows();row++){
2200 AliTPCclusterMI *cl = fOuterSec[sec][row].fClusters1;
2201 for (Int_t icl =0;icl< fOuterSec[sec][row].fN1;icl++)
2202 //if (cl[icl].IsUsed(10))
2204 cl = fOuterSec[sec][row].fClusters2;
2205 for (Int_t icl =0;icl< fOuterSec[sec][row].fN2;icl++)
2206 //if (cl[icl].IsUsed(10))
2215 void AliTPCtrackerMI::SignClusters(TObjArray * arr, Float_t fnumber, Float_t fdensity)
2218 //sign clusters to be "used"
2220 // snumber and sdensity sign number of sigmas - bellow mean value to be accepted
2221 // loop over "primaries"
2235 Int_t nseed = arr->GetEntriesFast();
2236 for (Int_t i=0; i<nseed; i++) {
2237 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2241 if (!(pt->IsActive())) continue;
2242 Float_t dens = pt->GetNumberOfClusters()/Float_t(pt->fNFoundable);
2243 if ( (dens>0.7) && (pt->GetNumberOfClusters()>70)){
2245 sumdens2+= dens*dens;
2246 sumn += pt->GetNumberOfClusters();
2247 sumn2 += pt->GetNumberOfClusters()*pt->GetNumberOfClusters();
2248 Float_t chi2 = pt->GetChi2()/pt->GetNumberOfClusters();
2251 sumchi2 +=chi2*chi2;
2256 Float_t mdensity = 0.9;
2257 Float_t meann = 130;
2258 Float_t meanchi = 1;
2259 Float_t sdensity = 0.1;
2260 Float_t smeann = 10;
2261 Float_t smeanchi =0.4;
2265 mdensity = sumdens/sum;
2267 meanchi = sumchi/sum;
2269 sdensity = sumdens2/sum-mdensity*mdensity;
2270 sdensity = TMath::Sqrt(sdensity);
2272 smeann = sumn2/sum-meann*meann;
2273 smeann = TMath::Sqrt(smeann);
2275 smeanchi = sumchi2/sum - meanchi*meanchi;
2276 smeanchi = TMath::Sqrt(smeanchi);
2280 //REMOVE SHORT DELTAS or tracks going out of sensitive volume of TPC
2282 for (Int_t i=0; i<nseed; i++) {
2283 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2287 if (pt->fBSigned) continue;
2288 if (pt->fBConstrain) continue;
2289 //if (!(pt->IsActive())) continue;
2291 Int_t found,foundable,shared;
2292 pt->GetClusterStatistic(0,160,found, foundable,shared);
2293 if (shared/float(found)>0.3) {
2294 if (shared/float(found)>0.9 ){
2295 //delete arr->RemoveAt(i);
2300 Bool_t isok =kFALSE;
2301 if ( (pt->fNShared/pt->GetNumberOfClusters()<0.5) &&pt->GetNumberOfClusters()>60)
2303 if ((TMath::Abs(1/pt->GetC())<100.) && (pt->fNShared/pt->GetNumberOfClusters()<0.7))
2305 if (TMath::Abs(pt->GetZ()/pt->GetX())>1.1)
2307 if ( (TMath::Abs(pt->GetSnp()>0.7) && pt->GetD(0,0)>60.))
2311 for (Int_t i=0; i<160; i++) {
2312 Int_t index=pt->GetClusterIndex2(i);
2313 if (index<0) continue;
2314 AliTPCclusterMI *c= pt->fClusterPointer[i];
2316 //if (!(c->IsUsed(10))) c->Use();
2323 Double_t maxchi = meanchi+2.*smeanchi;
2325 for (Int_t i=0; i<nseed; i++) {
2326 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2330 //if (!(pt->IsActive())) continue;
2331 if (pt->fBSigned) continue;
2332 Double_t chi = pt->GetChi2()/pt->GetNumberOfClusters();
2333 if (chi>maxchi) continue;
2336 Float_t dens = pt->GetNumberOfClusters()/Float_t(pt->fNFoundable);
2338 //sign only tracks with enoug big density at the beginning
2340 if ((pt->GetDensityFirst(40)<0.75) && pt->GetNumberOfClusters()<meann) continue;
2343 Double_t mindens = TMath::Max(double(mdensity-sdensity*fdensity*bfactor),0.65);
2344 Double_t minn = TMath::Max(Int_t(meann-fnumber*smeann*bfactor),50);
2346 // if (pt->fBConstrain) mindens = TMath::Max(mdensity-sdensity*fdensity*bfactor,0.65);
2347 if ( (pt->fRemoval==10) && (pt->GetSnp()>0.8)&&(dens>mindens))
2350 if ((dens>mindens && pt->GetNumberOfClusters()>minn) && chi<maxchi ){
2351 //Int_t noc=pt->GetNumberOfClusters();
2352 pt->fBSigned = kTRUE;
2353 for (Int_t i=0; i<160; i++) {
2355 Int_t index=pt->GetClusterIndex2(i);
2356 if (index<0) continue;
2357 AliTPCclusterMI *c= pt->fClusterPointer[i];
2359 // if (!(c->IsUsed(10))) c->Use();
2364 // gLastCheck = nseed;
2372 void AliTPCtrackerMI::StopNotActive(TObjArray * arr, Int_t row0, Float_t th0, Float_t th1, Float_t th2) const
2374 // stop not active tracks
2375 // take th1 as threshold for number of founded to number of foundable on last 10 active rows
2376 // take th2 as threshold for number of founded to number of foundable on last 20 active rows
2377 Int_t nseed = arr->GetEntriesFast();
2379 for (Int_t i=0; i<nseed; i++) {
2380 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2384 if (!(pt->IsActive())) continue;
2385 StopNotActive(pt,row0,th0, th1,th2);
2391 void AliTPCtrackerMI::StopNotActive(AliTPCseed * seed, Int_t row0, Float_t th0, Float_t th1,
2394 // stop not active tracks
2395 // take th1 as threshold for number of founded to number of foundable on last 10 active rows
2396 // take th2 as threshold for number of founded to number of foundable on last 20 active rows
2399 Int_t foundable = 0;
2400 Int_t maxindex = seed->fLastPoint; //last foundable row
2401 if (seed->fNFoundable*th0 > seed->GetNumberOfClusters()) {
2402 seed->Desactivate(10) ;
2406 for (Int_t i=row0; i<maxindex; i++){
2407 Int_t index = seed->GetClusterIndex2(i);
2408 if (index!=-1) foundable++;
2410 if (foundable<=30) sumgood1++;
2411 if (foundable<=50) {
2418 if (foundable>=30.){
2419 if (sumgood1<(th1*30.)) seed->Desactivate(10);
2422 if (sumgood2<(th2*50.)) seed->Desactivate(10);
2426 Int_t AliTPCtrackerMI::RefitInward(AliESD *event)
2429 // back propagation of ESD tracks
2435 //PrepareForProlongation(fSeeds,1);
2436 PropagateForward2(fSeeds);
2438 Int_t nseed = fSeeds->GetEntriesFast();
2439 for (Int_t i=0;i<nseed;i++){
2440 AliTPCseed * seed = (AliTPCseed*) fSeeds->UncheckedAt(i);
2441 if (!seed) continue;
2442 seed->PropagateTo(fParam->GetInnerRadiusLow());
2443 AliESDtrack *esd=event->GetTrack(i);
2444 seed->CookdEdx(0.02,0.6);
2445 CookLabel(seed,0.1); //For comparison only
2446 if (seed->GetNumberOfClusters()>60){
2447 esd->UpdateTrackParams(seed,AliESDtrack::kTPCrefit);
2451 //printf("problem\n");
2454 Info("RefitInward","Number of refitted tracks %d",ntracks);
2461 Int_t AliTPCtrackerMI::PropagateBack(AliESD *event)
2464 // back propagation of ESD tracks
2470 PropagateBack(fSeeds);
2471 Int_t nseed = fSeeds->GetEntriesFast();
2473 for (Int_t i=0;i<nseed;i++){
2474 AliTPCseed * seed = (AliTPCseed*) fSeeds->UncheckedAt(i);
2475 AliESDtrack *esd=event->GetTrack(i);
2476 seed->CookdEdx(0.02,0.6);
2477 CookLabel(seed,0.1); //For comparison only
2478 if (seed->GetNumberOfClusters()>60){
2479 esd->UpdateTrackParams(seed,AliESDtrack::kTPCout);
2483 Info("PropagateBack","Number of back propagated tracks %d",ntracks);
2490 void AliTPCtrackerMI::DeleteSeeds()
2494 Int_t nseed = fSeeds->GetEntriesFast();
2495 for (Int_t i=0;i<nseed;i++){
2496 AliTPCseed * seed = (AliTPCseed*)fSeeds->At(i);
2497 if (seed) delete fSeeds->RemoveAt(i);
2503 void AliTPCtrackerMI::ReadSeeds(AliESD *event, Int_t direction)
2506 //read seeds from the event
2508 Int_t nentr=event->GetNumberOfTracks();
2510 Info("ReadSeeds", "Number of ESD tracks: %d\n", nentr);
2515 fSeeds = new TObjArray(nentr);
2519 for (Int_t i=0; i<nentr; i++) {
2520 AliESDtrack *esd=event->GetTrack(i);
2521 ULong_t status=esd->GetStatus();
2522 AliTPCtrack t(*esd);
2523 AliTPCseed *seed = new AliTPCseed(t,t.GetAlpha());
2524 if ((status==AliESDtrack::kTPCin)&&(direction==1)) seed->ResetCovariance();
2525 if ( direction ==2 &&(status & AliESDtrack::kTRDrefit) == 0 ) seed->ResetCovariance();
2526 if ( direction ==2 && ((status & AliESDtrack::kTPCout) == 0) ) {
2531 if ( direction ==2 &&(status & AliESDtrack::kTRDrefit) > 0 ) {
2532 Double_t par0[5],par1[5],x;
2533 esd->GetInnerExternalParameters(x,par0);
2534 esd->GetExternalParameters(x,par1);
2535 Double_t delta1 = TMath::Abs(par0[4]-par1[4])/(0.000000001+TMath::Abs(par0[4]+par1[4]));
2536 Double_t delta2 = TMath::Abs(par0[3]-par1[3]);
2537 //reset covariance if suspicious
2538 if ( (delta1>0.1) || (delta2>0.006))
2539 seed->ResetCovariance();
2544 // rotate to the local coordinate system
2546 fSectors=fInnerSec; fN=fkNIS;
2548 Double_t alpha=seed->GetAlpha() - fSectors->GetAlphaShift();
2549 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
2550 if (alpha < 0. ) alpha += 2.*TMath::Pi();
2551 Int_t ns=Int_t(alpha/fSectors->GetAlpha())%fN;
2552 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
2553 if (alpha<-TMath::Pi()) alpha += 2*TMath::Pi();
2554 if (alpha>=TMath::Pi()) alpha -= 2*TMath::Pi();
2555 alpha-=seed->GetAlpha();
2556 if (!seed->Rotate(alpha)) {
2562 //seed->PropagateTo(fSectors->GetX(0));
2564 // Int_t index = esd->GetTPCindex();
2565 //AliTPCseed * seed2= (AliTPCseed*)fSeeds->At(index);
2566 //if (direction==2){
2567 // AliTPCseed * seed2 = ReSeed(seed,0.,0.5,1.);
2575 for (Int_t irow=0;irow<160;irow++){
2576 Int_t index = seed->GetClusterIndex2(irow);
2579 AliTPCclusterMI * cl = GetClusterMI(index);
2580 seed->fClusterPointer[irow] = cl;
2582 if ((index & 0x8000)==0){
2583 cl->Use(10); // accepted cluster
2585 cl->Use(6); // close cluster not accepted
2588 Info("ReadSeeds","Not found cluster");
2592 fSeeds->AddAt(seed,i);
2598 //_____________________________________________________________________________
2599 void AliTPCtrackerMI::MakeSeeds3(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4],
2600 Float_t deltay, Int_t ddsec) {
2601 //-----------------------------------------------------------------
2602 // This function creates track seeds.
2603 // SEEDING WITH VERTEX CONSTRAIN
2604 //-----------------------------------------------------------------
2605 // cuts[0] - fP4 cut
2606 // cuts[1] - tan(phi) cut
2607 // cuts[2] - zvertex cut
2608 // cuts[3] - fP3 cut
2616 Double_t x[5], c[15];
2617 // Int_t di = i1-i2;
2619 AliTPCseed * seed = new AliTPCseed;
2620 Double_t alpha=fSectors->GetAlpha(), shift=fSectors->GetAlphaShift();
2621 Double_t cs=cos(alpha), sn=sin(alpha);
2623 // Double_t x1 =fOuterSec->GetX(i1);
2624 //Double_t xx2=fOuterSec->GetX(i2);
2626 Double_t x1 =GetXrow(i1);
2627 Double_t xx2=GetXrow(i2);
2629 Double_t x3=GetX(), y3=GetY(), z3=GetZ();
2631 Int_t imiddle = (i2+i1)/2; //middle pad row index
2632 Double_t xm = GetXrow(imiddle); // radius of middle pad-row
2633 const AliTPCRow& krm=GetRow(sec,imiddle); //middle pad -row
2637 const AliTPCRow& kr1=GetRow(ns,i1);
2638 Double_t ymax = GetMaxY(i1)-kr1.fDeadZone-1.5;
2639 Double_t ymaxm = GetMaxY(imiddle)-kr1.fDeadZone-1.5;
2642 // change cut on curvature if it can't reach this layer
2643 // maximal curvature set to reach it
2644 Double_t dvertexmax = TMath::Sqrt((x1-x3)*(x1-x3)+(ymax+5-y3)*(ymax+5-y3));
2645 if (dvertexmax*0.5*cuts[0]>0.85){
2646 cuts[0] = 0.85/(dvertexmax*0.5+1.);
2648 Double_t r2min = 1/(cuts[0]*cuts[0]); //minimal square of radius given by cut
2651 if (deltay>0) ddsec = 0;
2652 // loop over clusters
2653 for (Int_t is=0; is < kr1; is++) {
2655 if (kr1[is]->IsUsed(10)) continue;
2656 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
2657 //if (TMath::Abs(y1)>ymax) continue;
2659 if (deltay>0 && TMath::Abs(ymax-TMath::Abs(y1))> deltay ) continue; // seed only at the edge
2661 // find possible directions
2662 Float_t anglez = (z1-z3)/(x1-x3);
2663 Float_t extraz = z1 - anglez*(x1-xx2); // extrapolated z
2666 //find rotation angles relative to line given by vertex and point 1
2667 Double_t dvertex2 = (x1-x3)*(x1-x3)+(y1-y3)*(y1-y3);
2668 Double_t dvertex = TMath::Sqrt(dvertex2);
2669 Double_t angle13 = TMath::ATan((y1-y3)/(x1-x3));
2670 Double_t cs13 = cos(-angle13), sn13 = sin(-angle13);
2673 // loop over 2 sectors
2679 Double_t dddz1=0; // direction of delta inclination in z axis
2686 for (Int_t dsec = dsec1; dsec<=dsec2;dsec++){
2687 Int_t sec2 = sec + dsec;
2689 // AliTPCRow& kr2 = fOuterSec[(sec2+fkNOS)%fkNOS][i2];
2690 //AliTPCRow& kr2m = fOuterSec[(sec2+fkNOS)%fkNOS][imiddle];
2691 AliTPCRow& kr2 = GetRow((sec2+fkNOS)%fkNOS,i2);
2692 AliTPCRow& kr2m = GetRow((sec2+fkNOS)%fkNOS,imiddle);
2693 Int_t index1 = TMath::Max(kr2.Find(extraz-0.6-dddz1*TMath::Abs(z1)*0.05)-1,0);
2694 Int_t index2 = TMath::Min(kr2.Find(extraz+0.6+dddz2*TMath::Abs(z1)*0.05)+1,kr2);
2696 // rotation angles to p1-p3
2697 Double_t cs13r = cos(-angle13+dsec*alpha)/dvertex, sn13r = sin(-angle13+dsec*alpha)/dvertex;
2698 Double_t x2, y2, z2;
2700 // Double_t dymax = maxangle*TMath::Abs(x1-xx2);
2703 Double_t dxx0 = (xx2-x3)*cs13r;
2704 Double_t dyy0 = (xx2-x3)*sn13r;
2705 for (Int_t js=index1; js < index2; js++) {
2706 const AliTPCclusterMI *kcl = kr2[js];
2707 if (kcl->IsUsed(10)) continue;
2709 //calcutate parameters
2711 Double_t yy0 = dyy0 +(kcl->GetY()-y3)*cs13r;
2713 if (TMath::Abs(yy0)<0.000001) continue;
2714 Double_t xx0 = dxx0 -(kcl->GetY()-y3)*sn13r;
2715 Double_t y0 = 0.5*(xx0*xx0+yy0*yy0-xx0)/yy0;
2716 Double_t r02 = (0.25+y0*y0)*dvertex2;
2717 //curvature (radius) cut
2718 if (r02<r2min) continue;
2722 Double_t c0 = 1/TMath::Sqrt(r02);
2726 //Double_t dfi0 = 2.*TMath::ASin(dvertex*c0*0.5);
2727 //Double_t dfi1 = 2.*TMath::ASin(TMath::Sqrt(yy0*yy0+(1-xx0)*(1-xx0))*dvertex*c0*0.5);
2728 Double_t dfi0 = 2.*AliTPCFastMath::FastAsin(dvertex*c0*0.5);
2729 Double_t dfi1 = 2.*AliTPCFastMath::FastAsin(TMath::Sqrt(yy0*yy0+(1-xx0)*(1-xx0))*dvertex*c0*0.5);
2732 Double_t z0 = kcl->GetZ();
2733 Double_t zzzz2 = z1-(z1-z3)*dfi1/dfi0;
2734 if (TMath::Abs(zzzz2-z0)>0.5) continue;
2737 Double_t dip = (z1-z0)*c0/dfi1;
2738 Double_t x0 = (0.5*cs13+y0*sn13)*dvertex*c0;
2749 x2= xx2*cs-y2*sn*dsec;
2750 y2=+xx2*sn*dsec+y2*cs;
2760 // do we have cluster at the middle ?
2762 GetProlongation(x1,xm,x,ym,zm);
2764 AliTPCclusterMI * cm=0;
2765 if (TMath::Abs(ym)-ymaxm<0){
2766 cm = krm.FindNearest2(ym,zm,1.0,0.6,dummy);
2767 if ((!cm) || (cm->IsUsed(10))) {
2772 // rotate y1 to system 0
2773 // get state vector in rotated system
2774 Double_t yr1 = (-0.5*sn13+y0*cs13)*dvertex*c0;
2775 Double_t xr2 = x0*cs+yr1*sn*dsec;
2776 Double_t xr[5]={kcl->GetY(),kcl->GetZ(), xr2, dip, c0};
2778 GetProlongation(xx2,xm,xr,ym,zm);
2779 if (TMath::Abs(ym)-ymaxm<0){
2780 cm = kr2m.FindNearest2(ym,zm,1.0,0.6,dummy);
2781 if ((!cm) || (cm->IsUsed(10))) {
2791 dym = ym - cm->GetY();
2792 dzm = zm - cm->GetZ();
2799 Double_t sy1=kr1[is]->GetSigmaY2()*2., sz1=kr1[is]->GetSigmaZ2()*2.;
2800 Double_t sy2=kcl->GetSigmaY2()*2., sz2=kcl->GetSigmaZ2()*2.;
2801 //Double_t sy3=400*3./12., sy=0.1, sz=0.1;
2802 Double_t sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
2803 //Double_t sy3=25000*x[4]*x[4]*60+0.5, sy=0.1, sz=0.1;
2805 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
2806 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
2807 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
2808 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
2809 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
2810 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
2812 Double_t f30=(F3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
2813 Double_t f31=(F3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
2814 Double_t f32=(F3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
2815 Double_t f34=(F3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
2819 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
2820 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
2821 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
2822 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
2823 c[13]=f30*sy1*f40+f32*sy2*f42;
2824 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
2826 // if (!BuildSeed(kr1[is],kcl,0,x1,x2,x3,x,c)) continue;
2828 UInt_t index=kr1.GetIndex(is);
2829 AliTPCseed *track=new(seed) AliTPCseed(index, x, c, x1, ns*alpha+shift);
2831 track->fIsSeeding = kTRUE;
2838 FollowProlongation(*track, (i1+i2)/2,1);
2839 Int_t foundable,found,shared;
2840 track->GetClusterStatistic((i1+i2)/2,i1, found, foundable, shared, kTRUE);
2841 if ((found<0.55*foundable) || shared>0.5*found || (track->GetSigmaY2()+track->GetSigmaZ2())>0.5){
2843 seed->~AliTPCseed();
2849 FollowProlongation(*track, i2,1);
2853 track->fBConstrain =1;
2854 // track->fLastPoint = i1+fInnerSec->GetNRows(); // first cluster in track position
2855 track->fLastPoint = i1; // first cluster in track position
2856 track->fFirstPoint = track->fLastPoint;
2858 if (track->GetNumberOfClusters()<(i1-i2)*0.5 ||
2859 track->GetNumberOfClusters() < track->fNFoundable*0.6 ||
2860 track->fNShared>0.4*track->GetNumberOfClusters() ) {
2862 seed->~AliTPCseed();
2866 // Z VERTEX CONDITION
2868 zv = track->GetZ()+track->GetTgl()/track->GetC()*
2869 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
2870 if (TMath::Abs(zv-z3)>cuts[2]) {
2871 FollowProlongation(*track, TMath::Max(i2-20,0));
2872 zv = track->GetZ()+track->GetTgl()/track->GetC()*
2873 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
2874 if (TMath::Abs(zv-z3)>cuts[2]){
2875 FollowProlongation(*track, TMath::Max(i2-40,0));
2876 zv = track->GetZ()+track->GetTgl()/track->GetC()*
2877 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
2878 if (TMath::Abs(zv-z3)>cuts[2] &&(track->GetNumberOfClusters() > track->fNFoundable*0.7)){
2879 // make seed without constrain
2880 AliTPCseed * track2 = MakeSeed(track,0.2,0.5,1.);
2881 FollowProlongation(*track2, i2,1);
2882 track2->fBConstrain = kFALSE;
2883 track2->fSeedType = 1;
2884 arr->AddLast(track2);
2886 seed->~AliTPCseed();
2891 seed->~AliTPCseed();
2898 track->fSeedType =0;
2899 arr->AddLast(track);
2900 seed = new AliTPCseed;
2902 // don't consider other combinations
2903 if (track->GetNumberOfClusters() > track->fNFoundable*0.8)
2909 Info("MakeSeeds3","\nSeeding statistic:\t%d\t%d\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin,nout1,nout2);
2915 void AliTPCtrackerMI::MakeSeeds5(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4],
2920 //-----------------------------------------------------------------
2921 // This function creates track seeds.
2922 //-----------------------------------------------------------------
2923 // cuts[0] - fP4 cut
2924 // cuts[1] - tan(phi) cut
2925 // cuts[2] - zvertex cut
2926 // cuts[3] - fP3 cut
2936 Double_t x[5], c[15];
2938 // make temporary seed
2939 AliTPCseed * seed = new AliTPCseed;
2940 Double_t alpha=fOuterSec->GetAlpha(), shift=fOuterSec->GetAlphaShift();
2941 // Double_t cs=cos(alpha), sn=sin(alpha);
2946 Double_t x1 = GetXrow(i1-1);
2947 const AliTPCRow& kr1=GetRow(sec,i1-1);
2948 Double_t y1max = GetMaxY(i1-1)-kr1.fDeadZone-1.5;
2950 Double_t x1p = GetXrow(i1);
2951 const AliTPCRow& kr1p=GetRow(sec,i1);
2953 Double_t x1m = GetXrow(i1-2);
2954 const AliTPCRow& kr1m=GetRow(sec,i1-2);
2957 //last 3 padrow for seeding
2958 AliTPCRow& kr3 = GetRow((sec+fkNOS)%fkNOS,i1-7);
2959 Double_t x3 = GetXrow(i1-7);
2960 // Double_t y3max= GetMaxY(i1-7)-kr3.fDeadZone-1.5;
2962 AliTPCRow& kr3p = GetRow((sec+fkNOS)%fkNOS,i1-6);
2963 Double_t x3p = GetXrow(i1-6);
2965 AliTPCRow& kr3m = GetRow((sec+fkNOS)%fkNOS,i1-8);
2966 Double_t x3m = GetXrow(i1-8);
2971 Int_t im = i1-4; //middle pad row index
2972 Double_t xm = GetXrow(im); // radius of middle pad-row
2973 const AliTPCRow& krm=GetRow(sec,im); //middle pad -row
2974 // Double_t ymmax = GetMaxY(im)-kr1.fDeadZone-1.5;
2977 Double_t deltax = x1-x3;
2978 Double_t dymax = deltax*cuts[1];
2979 Double_t dzmax = deltax*cuts[3];
2981 // loop over clusters
2982 for (Int_t is=0; is < kr1; is++) {
2984 if (kr1[is]->IsUsed(10)) continue;
2985 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
2987 if (deltay>0 && TMath::Abs(y1max-TMath::Abs(y1))> deltay ) continue; // seed only at the edge
2989 Int_t index1 = TMath::Max(kr3.Find(z1-dzmax)-1,0);
2990 Int_t index2 = TMath::Min(kr3.Find(z1+dzmax)+1,kr3);
2996 for (Int_t js=index1; js < index2; js++) {
2997 const AliTPCclusterMI *kcl = kr3[js];
2998 if (kcl->IsUsed(10)) continue;
3000 // apply angular cuts
3001 if (TMath::Abs(y1-y3)>dymax) continue;
3004 if (TMath::Abs(z1-z3)>dzmax) continue;
3006 Double_t angley = (y1-y3)/(x1-x3);
3007 Double_t anglez = (z1-z3)/(x1-x3);
3009 Double_t erry = TMath::Abs(angley)*(x1-x1m)*0.5+0.5;
3010 Double_t errz = TMath::Abs(anglez)*(x1-x1m)*0.5+0.5;
3012 Double_t yyym = angley*(xm-x1)+y1;
3013 Double_t zzzm = anglez*(xm-x1)+z1;
3015 const AliTPCclusterMI *kcm = krm.FindNearest2(yyym,zzzm,erry,errz,index);
3017 if (kcm->IsUsed(10)) continue;
3019 erry = TMath::Abs(angley)*(x1-x1m)*0.4+0.5;
3020 errz = TMath::Abs(anglez)*(x1-x1m)*0.4+0.5;
3027 // look around first
3028 const AliTPCclusterMI *kc1m = kr1m.FindNearest2(angley*(x1m-x1)+y1,
3034 if (kc1m->IsUsed(10)) used++;
3036 const AliTPCclusterMI *kc1p = kr1p.FindNearest2(angley*(x1p-x1)+y1,
3042 if (kc1p->IsUsed(10)) used++;
3044 if (used>1) continue;
3045 if (found<1) continue;
3049 const AliTPCclusterMI *kc3m = kr3m.FindNearest2(angley*(x3m-x3)+y3,
3055 if (kc3m->IsUsed(10)) used++;
3059 const AliTPCclusterMI *kc3p = kr3p.FindNearest2(angley*(x3p-x3)+y3,
3065 if (kc3p->IsUsed(10)) used++;
3069 if (used>1) continue;
3070 if (found<3) continue;
3080 x[4]=F1(x1,y1,x2,y2,x3,y3);
3081 //if (TMath::Abs(x[4]) >= cuts[0]) continue;
3084 x[2]=F2(x1,y1,x2,y2,x3,y3);
3087 x[3]=F3n(x1,y1,x2,y2,z1,z2,x[4]);
3088 //if (TMath::Abs(x[3]) > cuts[3]) continue;
3092 Double_t sy1=0.1, sz1=0.1;
3093 Double_t sy2=0.1, sz2=0.1;
3094 Double_t sy3=0.1, sy=0.1, sz=0.1;
3096 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
3097 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
3098 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
3099 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
3100 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
3101 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
3103 Double_t f30=(F3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
3104 Double_t f31=(F3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
3105 Double_t f32=(F3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
3106 Double_t f34=(F3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
3110 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3111 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3112 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3113 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3114 c[13]=f30*sy1*f40+f32*sy2*f42;
3115 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3117 // if (!BuildSeed(kr1[is],kcl,0,x1,x2,x3,x,c)) continue;
3119 UInt_t index=kr1.GetIndex(is);
3120 AliTPCseed *track=new(seed) AliTPCseed(index, x, c, x1, sec*alpha+shift);
3122 track->fIsSeeding = kTRUE;
3125 FollowProlongation(*track, i1-7,1);
3126 if (track->GetNumberOfClusters() < track->fNFoundable*0.75 ||
3127 track->fNShared>0.6*track->GetNumberOfClusters() || ( track->GetSigmaY2()+ track->GetSigmaZ2())>0.6){
3129 seed->~AliTPCseed();
3135 FollowProlongation(*track, i2,1);
3136 track->fBConstrain =0;
3137 track->fLastPoint = i1+fInnerSec->GetNRows(); // first cluster in track position
3138 track->fFirstPoint = track->fLastPoint;
3140 if (track->GetNumberOfClusters()<(i1-i2)*0.5 ||
3141 track->GetNumberOfClusters()<track->fNFoundable*0.7 ||
3142 track->fNShared>2. || track->GetChi2()/track->GetNumberOfClusters()>6 || ( track->GetSigmaY2()+ track->GetSigmaZ2())>0.5 ) {
3144 seed->~AliTPCseed();
3149 FollowProlongation(*track, TMath::Max(i2-10,0),1);
3150 AliTPCseed * track2 = MakeSeed(track,0.2,0.5,0.9);
3151 FollowProlongation(*track2, i2,1);
3152 track2->fBConstrain = kFALSE;
3153 track2->fSeedType = 4;
3154 arr->AddLast(track2);
3156 seed->~AliTPCseed();
3160 //arr->AddLast(track);
3161 //seed = new AliTPCseed;
3167 Info("MakeSeeds5","\nSeeding statiistic:\t%d\t%d\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin,nout1,nout2,nout3);
3173 //_____________________________________________________________________________
3174 void AliTPCtrackerMI::MakeSeeds2(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t */*cuts[4]*/,
3175 Float_t deltay, Bool_t /*bconstrain*/) {
3176 //-----------------------------------------------------------------
3177 // This function creates track seeds - without vertex constraint
3178 //-----------------------------------------------------------------
3179 // cuts[0] - fP4 cut - not applied
3180 // cuts[1] - tan(phi) cut
3181 // cuts[2] - zvertex cut - not applied
3182 // cuts[3] - fP3 cut
3192 Double_t alpha=fOuterSec->GetAlpha(), shift=fOuterSec->GetAlphaShift();
3193 // Double_t cs=cos(alpha), sn=sin(alpha);
3194 Int_t row0 = (i1+i2)/2;
3195 Int_t drow = (i1-i2)/2;
3196 const AliTPCRow& kr0=fSectors[sec][row0];
3199 AliTPCpolyTrack polytrack;
3200 Int_t nclusters=fSectors[sec][row0];
3201 AliTPCseed * seed = new AliTPCseed;
3206 for (Int_t is=0; is < nclusters; is++) { //LOOP over clusters
3208 Int_t nfoundable =0;
3209 for (Int_t iter =1; iter<2; iter++){ //iterations
3210 const AliTPCRow& krm=fSectors[sec][row0-iter];
3211 const AliTPCRow& krp=fSectors[sec][row0+iter];
3212 const AliTPCclusterMI * cl= kr0[is];
3214 if (cl->IsUsed(10)) {
3220 Double_t x = kr0.GetX();
3221 // Initialization of the polytrack
3226 Double_t y0= cl->GetY();
3227 Double_t z0= cl->GetZ();
3231 Double_t ymax = fSectors->GetMaxY(row0)-kr0.fDeadZone-1.5;
3232 if (deltay>0 && TMath::Abs(ymax-TMath::Abs(y0))> deltay ) continue; // seed only at the edge
3234 erry = (0.5)*cl->GetSigmaY2()/TMath::Sqrt(cl->GetQ())*6;
3235 errz = (0.5)*cl->GetSigmaZ2()/TMath::Sqrt(cl->GetQ())*6;
3236 polytrack.AddPoint(x,y0,z0,erry, errz);
3239 if (cl->IsUsed(10)) sumused++;
3242 Float_t roady = (5*TMath::Sqrt(cl->GetSigmaY2()+0.2)+1.)*iter;
3243 Float_t roadz = (5*TMath::Sqrt(cl->GetSigmaZ2()+0.2)+1.)*iter;
3246 AliTPCclusterMI * cl1 = krm.FindNearest(y0,z0,roady,roadz);
3247 if (cl1 && TMath::Abs(ymax-TMath::Abs(y0))) {
3248 erry = (0.5)*cl1->GetSigmaY2()/TMath::Sqrt(cl1->GetQ())*3;
3249 errz = (0.5)*cl1->GetSigmaZ2()/TMath::Sqrt(cl1->GetQ())*3;
3250 if (cl1->IsUsed(10)) sumused++;
3251 polytrack.AddPoint(x,cl1->GetY(),cl1->GetZ(),erry,errz);
3255 AliTPCclusterMI * cl2 = krp.FindNearest(y0,z0,roady,roadz);
3257 erry = (0.5)*cl2->GetSigmaY2()/TMath::Sqrt(cl2->GetQ())*3;
3258 errz = (0.5)*cl2->GetSigmaZ2()/TMath::Sqrt(cl2->GetQ())*3;
3259 if (cl2->IsUsed(10)) sumused++;
3260 polytrack.AddPoint(x,cl2->GetY(),cl2->GetZ(),erry,errz);
3263 if (sumused>0) continue;
3265 polytrack.UpdateParameters();
3271 nfoundable = polytrack.GetN();
3272 nfound = nfoundable;
3274 for (Int_t ddrow = iter+1; ddrow<drow;ddrow++){
3275 Float_t maxdist = 0.8*(1.+3./(ddrow));
3276 for (Int_t delta = -1;delta<=1;delta+=2){
3277 Int_t row = row0+ddrow*delta;
3278 kr = &(fSectors[sec][row]);
3279 Double_t xn = kr->GetX();
3280 Double_t ymax = fSectors->GetMaxY(row)-kr->fDeadZone-1.5;
3281 polytrack.GetFitPoint(xn,yn,zn);
3282 if (TMath::Abs(yn)>ymax) continue;
3284 AliTPCclusterMI * cln = kr->FindNearest(yn,zn,roady,roadz);
3286 Float_t dist = TMath::Sqrt( (yn-cln->GetY())*(yn-cln->GetY())+(zn-cln->GetZ())*(zn-cln->GetZ()));
3289 erry = (dist+0.3)*cln->GetSigmaY2()/TMath::Sqrt(cln->GetQ())*(1.+1./(ddrow));
3290 errz = (dist+0.3)*cln->GetSigmaZ2()/TMath::Sqrt(cln->GetQ())*(1.+1./(ddrow));
3291 if (cln->IsUsed(10)) {
3292 // printf("used\n");
3300 polytrack.AddPoint(xn,cln->GetY(),cln->GetZ(),erry, errz);
3305 if ( (sumused>3) || (sumused>0.5*nfound) || (nfound<0.6*nfoundable)) break;
3306 polytrack.UpdateParameters();
3309 if ( (sumused>3) || (sumused>0.5*nfound)) {
3310 //printf("sumused %d\n",sumused);
3315 polytrack.GetFitDerivation(kr0.GetX(),dy,dz);
3316 AliTPCpolyTrack track2;
3318 polytrack.Refit(track2,0.5+TMath::Abs(dy)*0.3,0.4+TMath::Abs(dz)*0.3);
3319 if (track2.GetN()<0.5*nfoundable) continue;
3322 if ((nfound>0.6*nfoundable) &&( nfoundable>0.4*(i1-i2))) {
3324 // test seed with and without constrain
3325 for (Int_t constrain=0; constrain<=0;constrain++){
3326 // add polytrack candidate
3328 Double_t x[5], c[15];
3329 Double_t x1,x2,x3,y1,y2,y3,z1,z2,z3;
3330 track2.GetBoundaries(x3,x1);
3332 track2.GetFitPoint(x1,y1,z1);
3333 track2.GetFitPoint(x2,y2,z2);
3334 track2.GetFitPoint(x3,y3,z3);
3336 //is track pointing to the vertex ?
3339 polytrack.GetFitPoint(x0,y0,z0);
3352 x[4]=F1(x1,y1,x2,y2,x3,y3);
3354 // if (TMath::Abs(x[4]) >= cuts[0]) continue; //
3355 x[2]=F2(x1,y1,x2,y2,x3,y3);
3357 //if (TMath::Abs(x[4]*x1-x[2]) >= cuts[1]) continue;
3358 //x[3]=F3(x1,y1,x2,y2,z1,z2);
3359 x[3]=F3n(x1,y1,x3,y3,z1,z3,x[4]);
3360 //if (TMath::Abs(x[3]) > cuts[3]) continue;
3363 Double_t sy =0.1, sz =0.1;
3364 Double_t sy1=0.02, sz1=0.02;
3365 Double_t sy2=0.02, sz2=0.02;
3369 sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
3372 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
3373 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
3374 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
3375 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
3376 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
3377 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
3379 Double_t f30=(F3(x1,y1+sy,x3,y3,z1,z3)-x[3])/sy;
3380 Double_t f31=(F3(x1,y1,x3,y3,z1+sz,z3)-x[3])/sz;
3381 Double_t f32=(F3(x1,y1,x3,y3+sy,z1,z3)-x[3])/sy;
3382 Double_t f34=(F3(x1,y1,x3,y3,z1,z3+sz)-x[3])/sz;
3387 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3388 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3389 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3390 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3391 c[13]=f30*sy1*f40+f32*sy2*f42;
3392 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3394 //Int_t row1 = fSectors->GetRowNumber(x1);
3395 Int_t row1 = GetRowNumber(x1);
3399 AliTPCseed *track=new (seed) AliTPCseed(index, x, c, x1, sec*alpha+shift);
3400 track->fIsSeeding = kTRUE;
3401 Int_t rc=FollowProlongation(*track, i2);
3402 if (constrain) track->fBConstrain =1;
3404 track->fBConstrain =0;
3405 track->fLastPoint = row1+fInnerSec->GetNRows(); // first cluster in track position
3406 track->fFirstPoint = track->fLastPoint;
3408 if (rc==0 || track->GetNumberOfClusters()<(i1-i2)*0.5 ||
3409 track->GetNumberOfClusters() < track->fNFoundable*0.6 ||
3410 track->fNShared>0.4*track->GetNumberOfClusters()) {
3413 seed->~AliTPCseed();
3416 arr->AddLast(track);
3417 seed = new AliTPCseed;
3421 } // if accepted seed
3424 Info("MakeSeeds2","\nSeeding statiistic:\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin3);
3430 AliTPCseed *AliTPCtrackerMI::MakeSeed(AliTPCseed *track, Float_t r0, Float_t r1, Float_t r2)
3434 //reseed using track points
3435 Int_t p0 = int(r0*track->GetNumberOfClusters()); // point 0
3436 Int_t p1 = int(r1*track->GetNumberOfClusters());
3437 Int_t p2 = int(r2*track->GetNumberOfClusters()); // last point
3439 Double_t x0[3],x1[3],x2[3];
3444 // find track position at given ratio of the length
3445 Int_t sec0, sec1, sec2;
3451 for (Int_t i=0;i<160;i++){
3452 if (track->fClusterPointer[i]){
3454 AliTPCTrackerPoint *trpoint =track->GetTrackPoint(i);
3455 if ( (index<p0) || x0[0]<0 ){
3456 if (trpoint->GetX()>1){
3457 clindex = track->GetClusterIndex2(i);
3459 x0[0] = trpoint->GetX();
3460 x0[1] = trpoint->GetY();
3461 x0[2] = trpoint->GetZ();
3462 sec0 = ((clindex&0xff000000)>>24)%18;
3467 if ( (index<p1) &&(trpoint->GetX()>1)){
3468 clindex = track->GetClusterIndex2(i);
3470 x1[0] = trpoint->GetX();
3471 x1[1] = trpoint->GetY();
3472 x1[2] = trpoint->GetZ();
3473 sec1 = ((clindex&0xff000000)>>24)%18;
3476 if ( (index<p2) &&(trpoint->GetX()>1)){
3477 clindex = track->GetClusterIndex2(i);
3479 x2[0] = trpoint->GetX();
3480 x2[1] = trpoint->GetY();
3481 x2[2] = trpoint->GetZ();
3482 sec2 = ((clindex&0xff000000)>>24)%18;
3489 Double_t alpha, cs,sn, xx2,yy2;
3491 alpha = (sec1-sec2)*fSectors->GetAlpha();
3492 cs = TMath::Cos(alpha);
3493 sn = TMath::Sin(alpha);
3494 xx2= x1[0]*cs-x1[1]*sn;
3495 yy2= x1[0]*sn+x1[1]*cs;
3499 alpha = (sec0-sec2)*fSectors->GetAlpha();
3500 cs = TMath::Cos(alpha);
3501 sn = TMath::Sin(alpha);
3502 xx2= x0[0]*cs-x0[1]*sn;
3503 yy2= x0[0]*sn+x0[1]*cs;
3509 Double_t x[5],c[15];
3513 x[4]=F1(x2[0],x2[1],x1[0],x1[1],x0[0],x0[1]);
3514 // if (x[4]>1) return 0;
3515 x[2]=F2(x2[0],x2[1],x1[0],x1[1],x0[0],x0[1]);
3516 x[3]=F3n(x2[0],x2[1],x0[0],x0[1],x2[2],x0[2],x[4]);
3517 //if (TMath::Abs(x[3]) > 2.2) return 0;
3518 //if (TMath::Abs(x[2]) > 1.99) return 0;
3520 Double_t sy =0.1, sz =0.1;
3522 Double_t sy1=0.02+track->GetSigmaY2(), sz1=0.02+track->GetSigmaZ2();
3523 Double_t sy2=0.01+track->GetSigmaY2(), sz2=0.01+track->GetSigmaZ2();
3524 Double_t sy3=0.01+track->GetSigmaY2();
3526 Double_t f40=(F1(x2[0],x2[1]+sy,x1[0],x1[1],x0[0],x0[1])-x[4])/sy;
3527 Double_t f42=(F1(x2[0],x2[1],x1[0],x1[1]+sy,x0[0],x0[1])-x[4])/sy;
3528 Double_t f43=(F1(x2[0],x2[1],x1[0],x1[1],x0[0],x0[1]+sy)-x[4])/sy;
3529 Double_t f20=(F2(x2[0],x2[1]+sy,x1[0],x1[1],x0[0],x0[1])-x[2])/sy;
3530 Double_t f22=(F2(x2[0],x2[1],x1[0],x1[1]+sy,x0[0],x0[1])-x[2])/sy;
3531 Double_t f23=(F2(x2[0],x2[1],x1[0],x1[1],x0[0],x0[1]+sy)-x[2])/sy;
3533 Double_t f30=(F3(x2[0],x2[1]+sy,x0[0],x0[1],x2[2],x0[2])-x[3])/sy;
3534 Double_t f31=(F3(x2[0],x2[1],x0[0],x0[1],x2[2]+sz,x0[2])-x[3])/sz;
3535 Double_t f32=(F3(x2[0],x2[1],x0[0],x0[1]+sy,x2[2],x0[2])-x[3])/sy;
3536 Double_t f34=(F3(x2[0],x2[1],x0[0],x0[1],x2[2],x0[2]+sz)-x[3])/sz;
3541 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3542 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3543 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3544 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3545 c[13]=f30*sy1*f40+f32*sy2*f42;
3546 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3548 // Int_t row1 = fSectors->GetRowNumber(x2[0]);
3549 AliTPCseed *seed=new AliTPCseed(0, x, c, x2[0], sec2*fSectors->GetAlpha()+fSectors->GetAlphaShift());
3550 // Double_t y0,z0,y1,z1, y2,z2;
3551 //seed->GetProlongation(x0[0],y0,z0);
3552 // seed->GetProlongation(x1[0],y1,z1);
3553 //seed->GetProlongation(x2[0],y2,z2);
3555 seed->fLastPoint = pp2;
3556 seed->fFirstPoint = pp2;
3563 AliTPCseed *AliTPCtrackerMI::ReSeed(AliTPCseed *track, Float_t r0, Float_t r1, Float_t r2)
3567 //reseed using founded clusters
3569 // Find the number of clusters
3570 Int_t nclusters = 0;
3571 for (Int_t irow=0;irow<160;irow++){
3572 if (track->GetClusterIndex(irow)>0) nclusters++;
3576 ipos[0] = TMath::Max(int(r0*nclusters),0); // point 0 cluster
3577 ipos[1] = TMath::Min(int(r1*nclusters),nclusters-1); //
3578 ipos[2] = TMath::Min(int(r2*nclusters),nclusters-1); // last point
3582 Int_t row[3],sec[3]={0,0,0};
3584 // find track row position at given ratio of the length
3586 for (Int_t irow=0;irow<160;irow++){
3587 if (track->GetClusterIndex2(irow)<0) continue;
3589 for (Int_t ipoint=0;ipoint<3;ipoint++){
3590 if (index<=ipos[ipoint]) row[ipoint] = irow;
3594 //Get cluster and sector position
3595 for (Int_t ipoint=0;ipoint<3;ipoint++){
3596 Int_t clindex = track->GetClusterIndex2(row[ipoint]);
3597 AliTPCclusterMI * cl = GetClusterMI(clindex);
3600 // AliTPCclusterMI * cl = GetClusterMI(clindex);
3603 sec[ipoint] = ((clindex&0xff000000)>>24)%18;
3604 xyz[ipoint][0] = GetXrow(row[ipoint]);
3605 xyz[ipoint][1] = cl->GetY();
3606 xyz[ipoint][2] = cl->GetZ();
3610 // Calculate seed state vector and covariance matrix
3612 Double_t alpha, cs,sn, xx2,yy2;
3614 alpha = (sec[1]-sec[2])*fSectors->GetAlpha();
3615 cs = TMath::Cos(alpha);
3616 sn = TMath::Sin(alpha);
3617 xx2= xyz[1][0]*cs-xyz[1][1]*sn;
3618 yy2= xyz[1][0]*sn+xyz[1][1]*cs;
3622 alpha = (sec[0]-sec[2])*fSectors->GetAlpha();
3623 cs = TMath::Cos(alpha);
3624 sn = TMath::Sin(alpha);
3625 xx2= xyz[0][0]*cs-xyz[0][1]*sn;
3626 yy2= xyz[0][0]*sn+xyz[0][1]*cs;
3632 Double_t x[5],c[15];
3636 x[4]=F1(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1]);
3637 x[2]=F2(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1]);
3638 x[3]=F3n(xyz[2][0],xyz[2][1],xyz[0][0],xyz[0][1],xyz[2][2],xyz[0][2],x[4]);
3640 Double_t sy =0.1, sz =0.1;
3642 Double_t sy1=0.2, sz1=0.2;
3643 Double_t sy2=0.2, sz2=0.2;
3646 Double_t f40=(F1(xyz[2][0],xyz[2][1]+sy,xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1])-x[4])/sy;
3647 Double_t f42=(F1(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1]+sy,xyz[0][0],xyz[0][1])-x[4])/sy;
3648 Double_t f43=(F1(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1]+sy)-x[4])/sy;
3649 Double_t f20=(F2(xyz[2][0],xyz[2][1]+sy,xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1])-x[2])/sy;
3650 Double_t f22=(F2(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1]+sy,xyz[0][0],xyz[0][1])-x[2])/sy;
3651 Double_t f23=(F2(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1]+sy)-x[2])/sy;
3653 Double_t f30=(F3(xyz[2][0],xyz[2][1]+sy,xyz[0][0],xyz[0][1],xyz[2][2],xyz[0][2])-x[3])/sy;
3654 Double_t f31=(F3(xyz[2][0],xyz[2][1],xyz[0][0],xyz[0][1],xyz[2][2]+sz,xyz[0][2])-x[3])/sz;
3655 Double_t f32=(F3(xyz[2][0],xyz[2][1],xyz[0][0],xyz[0][1]+sy,xyz[2][2],xyz[0][2])-x[3])/sy;
3656 Double_t f34=(F3(xyz[2][0],xyz[2][1],xyz[0][0],xyz[0][1],xyz[2][2],xyz[0][2]+sz)-x[3])/sz;
3661 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3662 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3663 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3664 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3665 c[13]=f30*sy1*f40+f32*sy2*f42;
3666 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3668 // Int_t row1 = fSectors->GetRowNumber(xyz[2][0]);
3669 AliTPCseed *seed=new AliTPCseed(0, x, c, xyz[2][0], sec[2]*fSectors->GetAlpha()+fSectors->GetAlphaShift());
3670 seed->fLastPoint = row[2];
3671 seed->fFirstPoint = row[2];
3675 Int_t AliTPCtrackerMI::CheckKinkPoint(AliTPCseed*seed, Float_t th)
3680 for (Int_t i=0;i<12;i++) seed->fKinkPoint[i]=0;
3682 if (TMath::Abs(seed->GetC())>0.01) return 0;
3685 Float_t x[160], y[160], erry[160], z[160], errz[160];
3687 Float_t xt[160], yt[160], zt[160];
3692 Int_t middle = seed->GetNumberOfClusters()/2;
3695 // find central sector, get local cooordinates
3697 for (Int_t i=seed->fFirstPoint;i<=seed->fLastPoint;i++) {
3698 sec[i]= seed->GetClusterSector(i)%18;
3701 AliTPCclusterMI * cl = seed->fClusterPointer[i];
3702 // if (cl==0) cl = GetClusterMI(seed->GetClusterIndex2(i));
3709 if (i>i2) i2 = i; //last point with cluster
3710 if (i2<i1) i1 = i; //first point with cluster
3713 AliTPCTrackerPoint * point = seed->GetTrackPoint(i);
3715 yt[i] = point->GetY();
3716 zt[i] = point->GetZ();
3718 if (point->GetX()>0){
3719 erry[i] = point->GetErrY();
3720 errz[i] = point->GetErrZ();
3725 secm = sec[i]; //central sector
3726 padm = i; //middle point with cluster
3731 // rotate position to global coordinate system connected to sector at last the point
3733 for (Int_t i=i1;i<=i2;i++){
3735 if (sec[i]<0) continue;
3736 Double_t alpha = (sec[i2]-sec[i])*fSectors->GetAlpha();
3737 Double_t cs = TMath::Cos(alpha);
3738 Double_t sn = TMath::Sin(alpha);
3739 Float_t xx2= x[i]*cs+y[i]*sn;
3740 Float_t yy2= -x[i]*sn+y[i]*cs;
3744 xx2= xt[i]*cs+yt[i]*sn;
3745 yy2= -xt[i]*sn+yt[i]*cs;
3750 //get "state" vector
3751 Double_t xh[5],xm = x[padm];
3754 xh[4]=F1(xt[i2],yt[i2],xt[padm],yt[padm],xt[i1],yt[i1]);
3755 xh[2]=F2(xt[i2],yt[i2],xt[padm],yt[padm],xt[i1],yt[i1]);
3756 xh[3]=F3n(xt[i2],yt[i2],xt[i1],yt[i1],zt[i2],zt[i1],xh[4]);
3759 for (Int_t i=i1;i<=i2;i++){
3761 if (sec[i]<0) continue;
3762 GetProlongation(x[i2], x[i],xh,yy,zz);
3763 if (TMath::Abs(y[i]-yy)>4||TMath::Abs(z[i]-zz)>4){
3765 //xxh[4]=F1old(x[i2],y[i2],x[padm],y[padm],x[i1],y[i1]);
3766 //xxh[2]=F2old(x[i2],y[i2],x[padm],y[padm],x[i1],y[i1]);
3767 Error("AliTPCtrackerMI::CheckKinkPoint","problem\n");
3772 Float_t dyup[160],dydown[160], dzup[160], dzdown[160];
3773 Float_t yup[160], ydown[160], zup[160], zdown[160];
3775 AliTPCpolyTrack ptrack1,ptrack2;
3778 for (Int_t i=i1;i<=i2;i++){
3779 AliTPCclusterMI * cl = seed->fClusterPointer[i];
3781 if (cl->GetType()<0) continue;
3782 if (cl->GetType()>10) continue;
3785 ptrack1.AddPoint(x[i]-xm,y[i],z[i],0.1,0.1);
3787 if (ptrack1.GetN()>4.){
3788 ptrack1.UpdateParameters();
3790 ptrack1.GetFitDerivation(x[i]-xm,ddy,ddz);
3792 ptrack1.GetFitPoint(x[i]-xm,yy,zz);
3801 dyup[i]=0.; //not enough points
3806 for (Int_t i=i2;i>=i1;i--){
3807 AliTPCclusterMI * cl = seed->fClusterPointer[i];
3809 if (cl->GetType()<0) continue;
3810 if (cl->GetType()>10) continue;
3812 ptrack2.AddPoint(x[i]-xm,y[i],z[i],0.1,0.1);
3814 if (ptrack2.GetN()>4){
3815 ptrack2.UpdateParameters();
3817 ptrack2.GetFitDerivation(x[i]-xm,ddy,ddz);
3819 ptrack2.GetFitPoint(x[i]-xm,yy,zz);
3827 dydown[i]=0.; //not enough points
3832 // find maximal difference of the derivation
3833 for (Int_t i=0;i<12;i++) seed->fKinkPoint[i]=0;
3836 for (Int_t i=i1+10;i<i2-10;i++){
3837 if ( (TMath::Abs(dydown[i])<0.00000001) || (TMath::Abs(dyup[i])<0.00000001) ||i<30)continue;
3838 // printf("%f\t%f\t%f\t%f\t%f\n",x[i],dydown[i],dyup[i],dzdown[i],dzup[i]);
3840 Float_t ddy = TMath::Abs(dydown[i]-dyup[i]);
3841 Float_t ddz = TMath::Abs(dzdown[i]-dzup[i]);
3842 if ( (ddy+ddz)> th){
3843 seed->fKinkPoint[0] = i;
3844 seed->fKinkPoint[1] = ddy;
3845 seed->fKinkPoint[2] = ddz;
3852 //write information to the debug tree
3853 TBranch * br = fTreeDebug->GetBranch("debug");
3854 TClonesArray * arr = new TClonesArray("AliTPCTrackPoint2");
3855 arr->ExpandCreateFast(i2-i1);
3856 br->SetAddress(&arr);
3858 AliTPCclusterMI cldummy;
3860 AliTPCTrackPoint2 pdummy;
3861 pdummy.GetTPoint().fIsShared = 10;
3863 Double_t alpha = sec[i2]*fSectors->GetAlpha();
3864 Double_t cs = TMath::Cos(alpha);
3865 Double_t sn = TMath::Sin(alpha);
3867 for (Int_t i=i1;i<i2;i++){
3868 AliTPCTrackPoint2 *trpoint = (AliTPCTrackPoint2*)arr->UncheckedAt(i-i1);
3870 AliTPCclusterMI * cl0 = seed->fClusterPointer[i];
3872 AliTPCTrackerPoint * point = seed->GetTrackPoint(i);
3875 Double_t x = GetXrow(i);
3876 trpoint->GetTPoint() = *point;
3877 trpoint->GetCPoint() = *cl0;
3878 trpoint->GetCPoint().SetQ(TMath::Abs(cl0->GetQ()));
3879 trpoint->fID = seed->GetUniqueID();
3880 trpoint->fLab = seed->GetLabel();
3882 trpoint->fGX = cs *x + sn*point->GetY();
3883 trpoint->fGY = -sn *x + cs*point->GetY() ;
3884 trpoint->fGZ = point->GetZ();
3886 trpoint->fDY = y[i];
3887 trpoint->fDZ = z[i];
3889 trpoint->fDYU = dyup[i];
3890 trpoint->fDZU = dzup[i];
3892 trpoint->fDYD = dydown[i];
3893 trpoint->fDZD = dzdown[i];
3895 if (TMath::Abs(dyup[i])>0.00000000001 &&TMath::Abs(dydown[i])>0.00000000001){
3896 trpoint->fDDY = dydown[i]-dyup[i];
3897 trpoint->fDDZ = dzdown[i]-dzup[i];
3905 trpoint->GetCPoint()= cldummy;
3922 AliTPCseed* AliTPCtrackerMI::ReSeed(AliTPCseed *t)
3925 // reseed - refit - track
3928 // Int_t last = fSectors->GetNRows()-1;
3930 if (fSectors == fOuterSec){
3931 first = TMath::Max(first, t->fFirstPoint-fInnerSec->GetNRows());
3935 first = t->fFirstPoint;
3937 AliTPCseed * seed = MakeSeed(t,0.1,0.5,0.9);
3938 FollowBackProlongation(*t,fSectors->GetNRows()-1);
3940 FollowProlongation(*t,first);
3950 //_____________________________________________________________________________
3951 Int_t AliTPCtrackerMI::ReadSeeds(const TFile *inp) {
3952 //-----------------------------------------------------------------
3953 // This function reades track seeds.
3954 //-----------------------------------------------------------------
3955 TDirectory *savedir=gDirectory;
3957 TFile *in=(TFile*)inp;
3958 if (!in->IsOpen()) {
3959 cerr<<"AliTPCtrackerMI::ReadSeeds(): input file is not open !\n";
3964 TTree *seedTree=(TTree*)in->Get("Seeds");
3966 cerr<<"AliTPCtrackerMI::ReadSeeds(): ";
3967 cerr<<"can't get a tree with track seeds !\n";
3970 AliTPCtrack *seed=new AliTPCtrack;
3971 seedTree->SetBranchAddress("tracks",&seed);
3973 if (fSeeds==0) fSeeds=new TObjArray(15000);
3975 Int_t n=(Int_t)seedTree->GetEntries();
3976 for (Int_t i=0; i<n; i++) {
3977 seedTree->GetEvent(i);
3978 fSeeds->AddLast(new AliTPCseed(*seed,seed->GetAlpha()));
3987 Int_t AliTPCtrackerMI::Clusters2Tracks (AliESD *esd)
3990 if (fSeeds) DeleteSeeds();
3993 if (!fSeeds) return 1;
4000 //_____________________________________________________________________________
4001 Int_t AliTPCtrackerMI::Clusters2Tracks() {
4002 //-----------------------------------------------------------------
4003 // This is a track finder.
4004 //-----------------------------------------------------------------
4005 TDirectory *savedir=gDirectory;
4009 fSeeds = Tracking();
4012 Info("Clusters2Tracks","Time for tracking: \t");timer.Print();timer.Start();
4014 //activate again some tracks
4015 for (Int_t i=0; i<fSeeds->GetEntriesFast(); i++) {
4016 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
4018 Int_t nc=t.GetNumberOfClusters();
4020 delete fSeeds->RemoveAt(i);
4023 if (pt->fRemoval==10) {
4024 if (pt->GetDensityFirst(20)>0.8 || pt->GetDensityFirst(30)>0.8 || pt->GetDensityFirst(40)>0.7)
4025 pt->Desactivate(10); // make track again active
4027 pt->Desactivate(20);
4028 delete fSeeds->RemoveAt(i);
4032 RemoveDouble(fSeeds,0.2,0.6,11);
4033 RemoveUsed(fSeeds,0.5,0.5,6);
4036 Int_t nseed=fSeeds->GetEntriesFast();
4038 for (Int_t i=0; i<nseed; i++) {
4039 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
4041 Int_t nc=t.GetNumberOfClusters();
4043 delete fSeeds->RemoveAt(i);
4046 CookLabel(pt,0.1); //For comparison only
4047 //if ((pt->IsActive() || (pt->fRemoval==10) )&& nc>50 &&pt->GetNumberOfClusters()>0.4*pt->fNFoundable){
4048 if ((pt->IsActive() || (pt->fRemoval==10) )){
4050 if (fDebug>0) cerr<<found<<'\r';
4054 delete fSeeds->RemoveAt(i);
4058 //RemoveOverlap(fSeeds,0.99,7,kTRUE);
4060 //RemoveUsed(fSeeds,0.9,0.9,6);
4062 nseed=fSeeds->GetEntriesFast();
4064 for (Int_t i=0; i<nseed; i++) {
4065 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
4067 Int_t nc=t.GetNumberOfClusters();
4069 delete fSeeds->RemoveAt(i);
4073 t.CookdEdx(0.02,0.6);
4074 // CheckKinkPoint(&t,0.05);
4075 //if ((pt->IsActive() || (pt->fRemoval==10) )&& nc>50 &&pt->GetNumberOfClusters()>0.4*pt->fNFoundable){
4076 if ((pt->IsActive() || (pt->fRemoval==10) )){
4084 delete fSeeds->RemoveAt(i);
4085 //AliTPCseed * seed1 = ReSeed(pt,0.05,0.5,1);
4087 // FollowProlongation(*seed1,0);
4088 // Int_t n = seed1->GetNumberOfClusters();
4089 // printf("fP4\t%f\t%f\n",seed1->GetC(),pt->GetC());
4090 // printf("fN\t%d\t%d\n", seed1->GetNumberOfClusters(),pt->GetNumberOfClusters());
4093 //AliTPCseed * seed2 = ReSeed(pt,0.95,0.5,0.05);
4097 SortTracks(fSeeds, 1);
4101 PrepareForBackProlongation(fSeeds,5.);
4102 PropagateBack(fSeeds);
4103 printf("Time for back propagation: \t");timer.Print();timer.Start();
4107 PrepareForProlongation(fSeeds,5.);
4108 PropagateForward2(fSeeds);
4110 printf("Time for FORWARD propagation: \t");timer.Print();timer.Start();
4111 // RemoveUsed(fSeeds,0.7,0.7,6);
4112 //RemoveOverlap(fSeeds,0.9,7,kTRUE);
4114 nseed=fSeeds->GetEntriesFast();
4116 for (Int_t i=0; i<nseed; i++) {
4117 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
4119 Int_t nc=t.GetNumberOfClusters();
4121 delete fSeeds->RemoveAt(i);
4124 t.CookdEdx(0.02,0.6);
4125 // CookLabel(pt,0.1); //For comparison only
4126 //if ((pt->IsActive() || (pt->fRemoval==10) )&& nc>50 &&pt->GetNumberOfClusters()>0.4*pt->fNFoundable){
4127 if ((pt->IsActive() || (pt->fRemoval==10) )){
4128 cerr<<found++<<'\r';
4131 delete fSeeds->RemoveAt(i);
4136 // fNTracks = found;
4138 Info("Clusters2Tracks","Time for overlap removal, track writing and dedx cooking: \t"); timer.Print();timer.Start();
4141 // cerr<<"Number of found tracks : "<<"\t"<<found<<endl;
4142 Info("Clusters2Tracks","Number of found tracks %d",found);
4144 // UnloadClusters();
4149 void AliTPCtrackerMI::Tracking(TObjArray * arr)
4152 // tracking of the seeds
4155 fSectors = fOuterSec;
4156 ParallelTracking(arr,150,63);
4157 fSectors = fOuterSec;
4158 ParallelTracking(arr,63,0);
4161 TObjArray * AliTPCtrackerMI::Tracking(Int_t seedtype, Int_t i1, Int_t i2, Float_t cuts[4], Float_t dy, Int_t dsec)
4166 TObjArray * arr = new TObjArray;
4168 fSectors = fOuterSec;
4171 for (Int_t sec=0;sec<fkNOS;sec++){
4172 if (seedtype==3) MakeSeeds3(arr,sec,i1,i2,cuts,dy, dsec);
4173 if (seedtype==4) MakeSeeds5(arr,sec,i1,i2,cuts,dy);
4174 if (seedtype==2) MakeSeeds2(arr,sec,i1,i2,cuts,dy);
4177 Info("Tracking","\nSeeding - %d\t%d\t%d\t%d\n",seedtype,i1,i2,arr->GetEntriesFast());
4189 TObjArray * AliTPCtrackerMI::Tracking()
4195 Int_t nup=fOuterSec->GetNRows()+fInnerSec->GetNRows();
4197 TObjArray * seeds = new TObjArray;
4206 Float_t fnumber = 3.0;
4207 Float_t fdensity = 3.0;
4212 for (Int_t delta = 0; delta<18; delta+=6){
4216 arr = Tracking(3,nup-1-delta,nup-1-delta-gap,cuts,-1,1);
4217 SumTracks(seeds,arr);
4218 SignClusters(seeds,fnumber,fdensity);
4220 for (Int_t i=2;i<6;i+=2){
4221 // seed high pt tracks
4224 arr = Tracking(3,nup-i-delta,nup-i-delta-gap,cuts,-1,0);
4225 SumTracks(seeds,arr);
4226 SignClusters(seeds,fnumber,fdensity);
4231 // RemoveUsed(seeds,0.9,0.9,1);
4232 // UnsignClusters();
4233 // SignClusters(seeds,fnumber,fdensity);
4237 for (Int_t delta = 20; delta<120; delta+=10){
4239 // seed high pt tracks
4243 arr = Tracking(3,nup-delta,nup-delta-gap,cuts,-1);
4244 SumTracks(seeds,arr);
4245 SignClusters(seeds,fnumber,fdensity);
4250 arr = Tracking(3,nup-delta-5,nup-delta-5-gap,cuts,-1);
4251 SumTracks(seeds,arr);
4252 SignClusters(seeds,fnumber,fdensity);
4263 Info("Tracking()","\n\nPrimary seeding\t%d\n\n",seeds->GetEntriesFast());
4267 // RemoveUsed(seeds,0.75,0.75,1);
4269 //SignClusters(seeds,fnumber,fdensity);
4278 arr = Tracking(4,nup-1,nup-1-gap,cuts,-1);
4279 SumTracks(seeds,arr);
4280 SignClusters(seeds,fnumber,fdensity);
4282 arr = Tracking(4,nup-2,nup-2-gap,cuts,-1);
4283 SumTracks(seeds,arr);
4284 SignClusters(seeds,fnumber,fdensity);
4286 arr = Tracking(4,nup-3,nup-3-gap,cuts,-1);
4287 SumTracks(seeds,arr);
4288 SignClusters(seeds,fnumber,fdensity);
4292 for (Int_t delta = 3; delta<30; delta+=5){
4298 arr = Tracking(4,nup-1-delta,nup-1-delta-gap,cuts,-1);
4299 SumTracks(seeds,arr);
4300 SignClusters(seeds,fnumber,fdensity);
4302 arr = Tracking(4,nup-3-delta,nup-5-delta-gap,cuts,4);
4303 SumTracks(seeds,arr);
4304 SignClusters(seeds,fnumber,fdensity);
4316 for (Int_t delta = 30; delta<70; delta+=10){
4322 arr = Tracking(4,nup-1-delta,nup-1-delta-gap,cuts,-1);
4323 SumTracks(seeds,arr);
4324 SignClusters(seeds,fnumber,fdensity);
4326 arr = Tracking(4,nup-5-delta,nup-5-delta-gap,cuts,5 );
4327 SumTracks(seeds,arr);
4328 SignClusters(seeds,fnumber,fdensity);
4332 Info("Tracking()","\n\nSecondary seeding\t%d\n\n",seeds->GetEntriesFast());
4343 void AliTPCtrackerMI::SumTracks(TObjArray *arr1,TObjArray *arr2) const
4346 //sum tracks to common container
4347 //remove suspicious tracks
4348 Int_t nseed = arr2->GetEntriesFast();
4349 for (Int_t i=0;i<nseed;i++){
4350 AliTPCseed *pt=(AliTPCseed*)arr2->UncheckedAt(i);
4353 // NORMAL ACTIVE TRACK
4354 if (pt->IsActive()){
4355 arr1->AddLast(arr2->RemoveAt(i));
4358 //remove not usable tracks
4359 if (pt->fRemoval!=10){
4360 delete arr2->RemoveAt(i);
4363 // REMOVE VERY SHORT TRACKS
4364 if (pt->GetNumberOfClusters()<20){
4365 delete arr2->RemoveAt(i);
4368 // ENABLE ONLY ENOUGH GOOD STOPPED TRACKS
4369 if (pt->GetDensityFirst(20)>0.8 || pt->GetDensityFirst(30)>0.8 || pt->GetDensityFirst(40)>0.7)
4370 arr1->AddLast(arr2->RemoveAt(i));
4372 delete arr2->RemoveAt(i);
4381 void AliTPCtrackerMI::ParallelTracking(TObjArray * arr, Int_t rfirst, Int_t rlast)
4384 // try to track in parralel
4386 Int_t nseed=arr->GetEntriesFast();
4387 //prepare seeds for tracking
4388 for (Int_t i=0; i<nseed; i++) {
4389 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i), &t=*pt;
4391 if (!t.IsActive()) continue;
4392 // follow prolongation to the first layer
4393 if ( (fSectors ==fInnerSec) || (t.fFirstPoint-fParam->GetNRowLow()>rfirst+1) )
4394 FollowProlongation(t, rfirst+1);
4399 for (Int_t nr=rfirst; nr>=rlast; nr--){
4400 if (nr<fInnerSec->GetNRows())
4401 fSectors = fInnerSec;
4403 fSectors = fOuterSec;
4404 // make indexes with the cluster tracks for given
4406 // find nearest cluster
4407 for (Int_t i=0; i<nseed; i++) {
4408 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i), &t=*pt;
4410 if (!pt->IsActive()) continue;
4411 // if ( (fSectors ==fOuterSec) && (pt->fFirstPoint-fParam->GetNRowLow())<nr) continue;
4412 if (pt->fRelativeSector>17) {
4415 UpdateClusters(t,nr);
4417 // prolonagate to the nearest cluster - if founded
4418 for (Int_t i=0; i<nseed; i++) {
4419 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
4421 if (!pt->IsActive()) continue;
4422 // if ((fSectors ==fOuterSec) && (pt->fFirstPoint-fParam->GetNRowLow())<nr) continue;
4423 if (pt->fRelativeSector>17) {
4426 FollowToNextCluster(*pt,nr);
4431 void AliTPCtrackerMI::PrepareForBackProlongation(TObjArray * arr,Float_t fac) const
4435 // if we use TPC track itself we have to "update" covariance
4437 Int_t nseed= arr->GetEntriesFast();
4438 for (Int_t i=0;i<nseed;i++){
4439 AliTPCseed *pt = (AliTPCseed*)arr->UncheckedAt(i);
4443 //rotate to current local system at first accepted point
4444 Int_t index = pt->GetClusterIndex2(pt->fFirstPoint);
4445 Int_t sec = (index&0xff000000)>>24;
4447 Float_t angle1 = fInnerSec->GetAlpha()*sec+fInnerSec->GetAlphaShift();
4448 if (angle1>TMath::Pi())
4449 angle1-=2.*TMath::Pi();
4450 Float_t angle2 = pt->GetAlpha();
4452 if (TMath::Abs(angle1-angle2)>0.001){
4453 pt->Rotate(angle1-angle2);
4454 //angle2 = pt->GetAlpha();
4455 //pt->fRelativeSector = pt->GetAlpha()/fInnerSec->GetAlpha();
4456 //if (pt->GetAlpha()<0)
4457 // pt->fRelativeSector+=18;
4458 //sec = pt->fRelativeSector;
4467 void AliTPCtrackerMI::PrepareForProlongation(TObjArray * arr, Float_t fac) const
4471 // if we use TPC track itself we have to "update" covariance
4473 Int_t nseed= arr->GetEntriesFast();
4474 for (Int_t i=0;i<nseed;i++){
4475 AliTPCseed *pt = (AliTPCseed*)arr->UncheckedAt(i);
4478 pt->fFirstPoint = pt->fLastPoint;
4486 Int_t AliTPCtrackerMI::PropagateBack(TObjArray * arr)
4489 // make back propagation
4491 Int_t nseed= arr->GetEntriesFast();
4492 for (Int_t i=0;i<nseed;i++){
4493 AliTPCseed *pt = (AliTPCseed*)arr->UncheckedAt(i);
4495 //AliTPCseed *pt2 = new AliTPCseed(*pt);
4496 fSectors = fInnerSec;
4497 //FollowBackProlongation(*pt,fInnerSec->GetNRows()-1);
4498 //fSectors = fOuterSec;
4499 FollowBackProlongation(*pt,fInnerSec->GetNRows()+fOuterSec->GetNRows()-1);
4500 //if (pt->GetNumberOfClusters()<(pt->fEsd->GetTPCclusters(0)) ){
4501 // Error("PropagateBack","Not prolonged track %d",pt->GetLabel());
4502 // FollowBackProlongation(*pt2,fInnerSec->GetNRows()+fOuterSec->GetNRows()-1);
4510 Int_t AliTPCtrackerMI::PropagateForward2(TObjArray * arr)
4513 // make forward propagation
4515 Int_t nseed= arr->GetEntriesFast();
4517 for (Int_t i=0;i<nseed;i++){
4518 AliTPCseed *pt = (AliTPCseed*)arr->UncheckedAt(i);
4520 FollowProlongation(*pt,0);
4527 Int_t AliTPCtrackerMI::PropagateForward()
4530 // propagate track forward
4532 Int_t nseed = fSeeds->GetEntriesFast();
4533 for (Int_t i=0;i<nseed;i++){
4534 AliTPCseed *pt = (AliTPCseed*)fSeeds->UncheckedAt(i);
4536 AliTPCseed &t = *pt;
4537 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
4538 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
4539 if (alpha < 0. ) alpha += 2.*TMath::Pi();
4540 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
4544 fSectors = fOuterSec;
4545 ParallelTracking(fSeeds,fOuterSec->GetNRows()+fInnerSec->GetNRows()-1,fInnerSec->GetNRows());
4546 fSectors = fInnerSec;
4547 ParallelTracking(fSeeds,fInnerSec->GetNRows()-1,0);
4557 Int_t AliTPCtrackerMI::PropagateBack(AliTPCseed * pt, Int_t row0, Int_t row1)
4560 // make back propagation, in between row0 and row1
4564 fSectors = fInnerSec;
4567 if (row1<fSectors->GetNRows())
4570 r1 = fSectors->GetNRows()-1;
4572 if (row0<fSectors->GetNRows()&& r1>0 )
4573 FollowBackProlongation(*pt,r1);
4574 if (row1<=fSectors->GetNRows())
4577 r1 = row1 - fSectors->GetNRows();
4578 if (r1<=0) return 0;
4579 if (r1>=fOuterSec->GetNRows()) return 0;
4580 fSectors = fOuterSec;
4581 return FollowBackProlongation(*pt,r1);
4589 void AliTPCtrackerMI::GetShape(AliTPCseed * seed, Int_t row)
4593 Float_t sd2 = TMath::Abs((fParam->GetZLength()-TMath::Abs(seed->GetZ())))*fParam->GetDiffL()*fParam->GetDiffL();
4594 // Float_t padlength = fParam->GetPadPitchLength(seed->fSector);
4595 Float_t padlength = GetPadPitchLength(row);
4597 Float_t sresy = (seed->fSector < fParam->GetNSector()/2) ? 0.2 :0.3;
4598 Float_t angulary = seed->GetSnp();
4599 angulary = angulary*angulary/(1-angulary*angulary);
4600 seed->fCurrentSigmaY2 = sd2+padlength*padlength*angulary/12.+sresy*sresy;
4602 Float_t sresz = fParam->GetZSigma();
4603 Float_t angularz = seed->GetTgl();
4604 seed->fCurrentSigmaZ2 = sd2+padlength*padlength*angularz*angularz*(1+angulary)/12.+sresz*sresz;
4606 Float_t wy = GetSigmaY(seed);
4607 Float_t wz = GetSigmaZ(seed);
4610 if (TMath::Abs(wy/seed->fCurrentSigmaY2-1)>0.0001 || TMath::Abs(wz/seed->fCurrentSigmaZ2-1)>0.0001 ){
4611 printf("problem\n");
4617 Float_t AliTPCtrackerMI::GetSigmaY(AliTPCseed * seed)
4621 Float_t sd2 = TMath::Abs((fParam->GetZLength()-TMath::Abs(seed->GetZ())))*fParam->GetDiffL()*fParam->GetDiffL();
4622 Float_t padlength = fParam->GetPadPitchLength(seed->fSector);
4623 Float_t sres = (seed->fSector < fParam->GetNSector()/2) ? 0.2 :0.3;
4624 Float_t angular = seed->GetSnp();
4625 angular = angular*angular/(1-angular*angular);
4626 // angular*=angular;
4627 //angular = TMath::Sqrt(angular/(1-angular));
4628 Float_t res = TMath::Sqrt(sd2+padlength*padlength*angular/12.+sres*sres);
4631 Float_t AliTPCtrackerMI::GetSigmaZ(AliTPCseed * seed)
4635 Float_t sd2 = TMath::Abs((fParam->GetZLength()-TMath::Abs(seed->GetZ())))*fParam->GetDiffL()*fParam->GetDiffL();
4636 Float_t padlength = fParam->GetPadPitchLength(seed->fSector);
4637 Float_t sres = fParam->GetZSigma();
4638 Float_t angular = seed->GetTgl();
4639 Float_t res = TMath::Sqrt(sd2+padlength*padlength*angular*angular/12.+sres*sres);
4646 //__________________________________________________________________________
4647 void AliTPCtrackerMI::CookLabel(AliTPCseed *t, Float_t wrong) const {
4648 //--------------------------------------------------------------------
4649 //This function "cooks" a track label. If label<0, this track is fake.
4650 //--------------------------------------------------------------------
4651 Int_t noc=t->GetNumberOfClusters();
4653 //printf("\nnot founded prolongation\n\n\n");
4659 AliTPCclusterMI *clusters[160];
4661 for (Int_t i=0;i<160;i++) {
4668 for (i=0; i<160 && current<noc; i++) {
4670 Int_t index=t->GetClusterIndex2(i);
4671 if (index<=0) continue;
4672 if (index&0x8000) continue;
4674 //clusters[current]=GetClusterMI(index);
4675 if (t->fClusterPointer[i]){
4676 clusters[current]=t->fClusterPointer[i];
4682 Int_t lab=123456789;
4683 for (i=0; i<noc; i++) {
4684 AliTPCclusterMI *c=clusters[i];
4686 lab=TMath::Abs(c->GetLabel(0));
4688 for (j=0; j<noc; j++) if (lb[j]==lab || mx[j]==0) break;
4694 for (i=0; i<noc; i++) if (mx[i]>max) {max=mx[i]; lab=lb[i];}
4696 for (i=0; i<noc; i++) {
4697 AliTPCclusterMI *c=clusters[i];
4699 if (TMath::Abs(c->GetLabel(1)) == lab ||
4700 TMath::Abs(c->GetLabel(2)) == lab ) max++;
4703 if ((1.- Float_t(max)/noc) > wrong) lab=-lab;
4706 Int_t tail=Int_t(0.10*noc);
4709 for (i=1; i<=160&&ind<tail; i++) {
4710 // AliTPCclusterMI *c=clusters[noc-i];
4711 AliTPCclusterMI *c=clusters[i];
4713 if (lab == TMath::Abs(c->GetLabel(0)) ||
4714 lab == TMath::Abs(c->GetLabel(1)) ||
4715 lab == TMath::Abs(c->GetLabel(2))) max++;
4718 if (max < Int_t(0.5*tail)) lab=-lab;
4725 //delete[] clusters;
4729 Int_t AliTPCtrackerMI::AliTPCSector::GetRowNumber(Double_t x) const
4731 //return pad row number for this x
4734 r=fRow[fN-1].GetX();
4735 if (x > r) return fN;
4737 if (x < r) return -1;
4738 return Int_t((x-r)/fPadPitchLength + 0.5);}
4740 r=fRow[fN-1].GetX();
4741 if (x > r) return fN;
4743 if (x < r) return -1;
4744 Double_t r1=fRow[64].GetX();
4746 return Int_t((x-r)/f1PadPitchLength + 0.5);}
4748 return (Int_t((x-r1)/f2PadPitchLength + 0.5)+64);}
4752 //_________________________________________________________________________
4753 void AliTPCtrackerMI::AliTPCSector::Setup(const AliTPCParam *par, Int_t f) {
4754 //-----------------------------------------------------------------------
4755 // Setup inner sector
4756 //-----------------------------------------------------------------------
4758 fAlpha=par->GetInnerAngle();
4759 fAlphaShift=par->GetInnerAngleShift();
4760 fPadPitchWidth=par->GetInnerPadPitchWidth();
4761 fPadPitchLength=par->GetInnerPadPitchLength();
4762 fN=par->GetNRowLow();
4763 fRow=new AliTPCRow[fN];
4764 for (Int_t i=0; i<fN; i++) {
4765 fRow[i].SetX(par->GetPadRowRadiiLow(i));
4766 fRow[i].fDeadZone =1.5; //1.5 cm of dead zone
4769 fAlpha=par->GetOuterAngle();
4770 fAlphaShift=par->GetOuterAngleShift();
4771 fPadPitchWidth = par->GetOuterPadPitchWidth();
4772 fPadPitchLength = par->GetOuter1PadPitchLength();
4773 f1PadPitchLength = par->GetOuter1PadPitchLength();
4774 f2PadPitchLength = par->GetOuter2PadPitchLength();
4776 fN=par->GetNRowUp();
4777 fRow=new AliTPCRow[fN];
4778 for (Int_t i=0; i<fN; i++) {
4779 fRow[i].SetX(par->GetPadRowRadiiUp(i));
4780 fRow[i].fDeadZone =1.5; // 1.5 cm of dead zone
4785 AliTPCtrackerMI::AliTPCRow::AliTPCRow() {
4787 // default constructor
4795 AliTPCtrackerMI::AliTPCRow::~AliTPCRow(){
4802 //_________________________________________________________________________
4804 AliTPCtrackerMI::AliTPCRow::InsertCluster(const AliTPCclusterMI* c, UInt_t index) {
4805 //-----------------------------------------------------------------------
4806 // Insert a cluster into this pad row in accordence with its y-coordinate
4807 //-----------------------------------------------------------------------
4808 if (fN==kMaxClusterPerRow) {
4809 cerr<<"AliTPCRow::InsertCluster(): Too many clusters !\n"; return;
4811 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
4812 Int_t i=Find(c->GetZ());
4813 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTPCclusterMI*));
4814 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
4815 fIndex[i]=index; fClusters[i]=c; fN++;
4818 void AliTPCtrackerMI::AliTPCRow::ResetClusters() {
4824 //delete[] fClusterArray;
4825 if (fClusters1) delete []fClusters1;
4826 if (fClusters2) delete []fClusters2;
4833 //___________________________________________________________________
4834 Int_t AliTPCtrackerMI::AliTPCRow::Find(Double_t z) const {
4835 //-----------------------------------------------------------------------
4836 // Return the index of the nearest cluster
4837 //-----------------------------------------------------------------------
4838 if (fN==0) return 0;
4839 if (z <= fClusters[0]->GetZ()) return 0;
4840 if (z > fClusters[fN-1]->GetZ()) return fN;
4841 Int_t b=0, e=fN-1, m=(b+e)/2;
4842 for (; b<e; m=(b+e)/2) {
4843 if (z > fClusters[m]->GetZ()) b=m+1;
4851 //___________________________________________________________________
4852 AliTPCclusterMI * AliTPCtrackerMI::AliTPCRow::FindNearest(Double_t y, Double_t z, Double_t roady, Double_t roadz) const {
4853 //-----------------------------------------------------------------------
4854 // Return the index of the nearest cluster in z y
4855 //-----------------------------------------------------------------------
4856 Float_t maxdistance = roady*roady + roadz*roadz;
4858 AliTPCclusterMI *cl =0;
4859 for (Int_t i=Find(z-roadz); i<fN; i++) {
4860 AliTPCclusterMI *c=(AliTPCclusterMI*)(fClusters[i]);
4861 if (c->GetZ() > z+roadz) break;
4862 if ( (c->GetY()-y) > roady ) continue;
4863 Float_t distance = (c->GetZ()-z)*(c->GetZ()-z)+(c->GetY()-y)*(c->GetY()-y);
4864 if (maxdistance>distance) {
4865 maxdistance = distance;
4872 AliTPCclusterMI * AliTPCtrackerMI::AliTPCRow::FindNearest2(Double_t y, Double_t z, Double_t roady, Double_t roadz,UInt_t & index) const
4874 //-----------------------------------------------------------------------
4875 // Return the index of the nearest cluster in z y
4876 //-----------------------------------------------------------------------
4877 Float_t maxdistance = roady*roady + roadz*roadz;
4878 Int_t iz1 = TMath::Max(fFastCluster[Int_t(z-roadz+254.5)]-1,0);
4879 Int_t iz2 = TMath::Min(fFastCluster[Int_t(z+roadz+255.5)]+1,fN);
4881 AliTPCclusterMI *cl =0;
4882 //FindNearest3(y,z,roady,roadz,index);
4883 // for (Int_t i=Find(z-roadz); i<fN; i++) {
4884 for (Int_t i=iz1; i<iz2; i++) {
4885 AliTPCclusterMI *c=(AliTPCclusterMI*)(fClusters[i]);
4886 if (c->GetZ() > z+roadz) break;
4887 if ( c->GetY()-y > roady ) continue;
4888 if ( y-c->GetY() > roady ) continue;
4889 Float_t distance = (c->GetZ()-z)*(c->GetZ()-z)+(c->GetY()-y)*(c->GetY()-y);
4890 if (maxdistance>distance) {
4891 maxdistance = distance;
4894 //roady = TMath::Sqrt(maxdistance);
4902 AliTPCclusterMI * AliTPCtrackerMI::AliTPCRow::FindNearest3(Double_t y, Double_t z, Double_t roady, Double_t roadz,UInt_t & index) const
4904 //-----------------------------------------------------------------------
4905 // Return the index of the nearest cluster in z y
4906 //-----------------------------------------------------------------------
4907 Float_t maxdistance = roady*roady + roadz*roadz;
4908 // Int_t iz = Int_t(z+255.);
4909 AliTPCclusterMI *cl =0;
4910 for (Int_t i=Find(z-roadz); i<fN; i++) {
4911 //for (Int_t i=fFastCluster[iz-2]; i<fFastCluster[iz+2]; i++) {
4912 AliTPCclusterMI *c=(AliTPCclusterMI*)(fClusters[i]);
4913 if (c->GetZ() > z+roadz) break;
4914 if ( c->GetY()-y > roady ) continue;
4915 if ( y-c->GetY() > roady ) continue;
4916 Float_t distance = (c->GetZ()-z)*(c->GetZ()-z)+(c->GetY()-y)*(c->GetY()-y);
4917 if (maxdistance>distance) {
4918 maxdistance = distance;
4921 //roady = TMath::Sqrt(maxdistance);
4930 AliTPCseed::AliTPCseed():AliTPCtrack(){
4934 for (Int_t i=0;i<200;i++) SetClusterIndex2(i,-3);
4935 for (Int_t i=0;i<160;i++) fClusterPointer[i]=0;
4953 AliTPCseed::AliTPCseed(const AliTPCtrack &t):AliTPCtrack(t){
4962 for (Int_t i=0;i<160;i++) {
4963 fClusterPointer[i] = 0;
4964 Int_t index = t.GetClusterIndex(i);
4966 SetClusterIndex2(i,index);
4969 SetClusterIndex2(i,-3);
4982 AliTPCseed::AliTPCseed(const AliKalmanTrack &t, Double_t a):AliTPCtrack(t,a){
4986 for (Int_t i=0;i<160;i++) {
4987 fClusterPointer[i] = 0;
4988 Int_t index = t.GetClusterIndex(i);
4989 SetClusterIndex2(i,index);
5010 AliTPCseed::AliTPCseed(UInt_t index, const Double_t xx[5], const Double_t cc[15],
5011 Double_t xr, Double_t alpha):
5012 AliTPCtrack(index, xx, cc, xr, alpha) {
5017 for (Int_t i=0;i<200;i++) SetClusterIndex2(i,-3);
5018 for (Int_t i=0;i<160;i++) fClusterPointer[i]=0;
5027 // fHelixIn = new TClonesArray("AliHelix",0);
5028 //fHelixOut = new TClonesArray("AliHelix",0);
5038 AliTPCseed::~AliTPCseed(){
5041 if (fPoints) delete fPoints;
5043 if (fEPoints) delete fEPoints;
5048 AliTPCTrackerPoint * AliTPCseed::GetTrackPoint(Int_t i)
5052 return &fTrackPoints[i];
5055 void AliTPCseed::RebuildSeed()
5058 // rebuild seed to be ready for storing
5059 AliTPCclusterMI cldummy;
5061 AliTPCTrackPoint pdummy;
5062 pdummy.GetTPoint().fIsShared = 10;
5063 for (Int_t i=0;i<160;i++){
5064 AliTPCclusterMI * cl0 = fClusterPointer[i];
5065 AliTPCTrackPoint *trpoint = (AliTPCTrackPoint*)fPoints->UncheckedAt(i);
5067 trpoint->GetTPoint() = *(GetTrackPoint(i));
5068 trpoint->GetCPoint() = *cl0;
5069 trpoint->GetCPoint().SetQ(TMath::Abs(cl0->GetQ()));
5073 trpoint->GetCPoint()= cldummy;
5081 Double_t AliTPCseed::GetDensityFirst(Int_t n)
5085 // return cluster for n rows bellow first point
5086 Int_t nfoundable = 1;
5088 for (Int_t i=fLastPoint-1;i>0&&nfoundable<n; i--){
5089 Int_t index = GetClusterIndex2(i);
5090 if (index!=-1) nfoundable++;
5091 if (index>0) nfound++;
5093 if (nfoundable<n) return 0;
5094 return Double_t(nfound)/Double_t(nfoundable);
5099 void AliTPCseed::GetClusterStatistic(Int_t first, Int_t last, Int_t &found, Int_t &foundable, Int_t &shared, Bool_t plus2)
5101 // get cluster stat. on given region
5106 for (Int_t i=first;i<last; i++){
5107 Int_t index = GetClusterIndex2(i);
5108 if (index!=-1) foundable++;
5109 if (fClusterPointer[i]) {
5115 if (fClusterPointer[i]->IsUsed(10)) {
5119 if (!plus2) continue; //take also neighborhoud
5121 if ( (i>0) && fClusterPointer[i-1]){
5122 if (fClusterPointer[i-1]->IsUsed(10)) {
5127 if ( fClusterPointer[i+1]){
5128 if (fClusterPointer[i+1]->IsUsed(10)) {
5135 //if (shared>found){
5136 //Error("AliTPCseed::GetClusterStatistic","problem\n");
5140 //_____________________________________________________________________________
5141 void AliTPCseed::CookdEdx(Double_t low, Double_t up,Int_t i1, Int_t i2, Bool_t onlyused) {
5142 //-----------------------------------------------------------------
5143 // This funtion calculates dE/dX within the "low" and "up" cuts.
5144 //-----------------------------------------------------------------
5147 Float_t angular[200];
5148 Float_t weight[200];
5151 // TClonesArray & arr = *fPoints;
5152 Float_t meanlog = 100.;
5154 Float_t mean[4] = {0,0,0,0};
5155 Float_t sigma[4] = {1000,1000,1000,1000};
5156 Int_t nc[4] = {0,0,0,0};
5157 Float_t norm[4] = {1000,1000,1000,1000};
5162 for (Int_t of =0; of<4; of++){
5163 for (Int_t i=of+i1;i<i2;i+=4)
5165 Int_t index = fIndex[i];
5166 if (index<0||index&0x8000) continue;
5168 //AliTPCTrackPoint * point = (AliTPCTrackPoint *) arr.At(i);
5169 AliTPCTrackerPoint * point = GetTrackPoint(i);
5170 //AliTPCTrackerPoint * pointm = GetTrackPoint(i-1);
5171 //AliTPCTrackerPoint * pointp = 0;
5172 //if (i<159) pointp = GetTrackPoint(i+1);
5174 if (point==0) continue;
5175 AliTPCclusterMI * cl = fClusterPointer[i];
5176 if (cl==0) continue;
5177 if (onlyused && (!cl->IsUsed(10))) continue;
5178 if (cl->IsUsed(11)) {
5182 Int_t type = cl->GetType();
5183 //if (point->fIsShared){
5188 // if (pointm->fIsShared) continue;
5190 // if (pointp->fIsShared) continue;
5192 if (type<0) continue;
5193 //if (type>10) continue;
5194 //if (point->GetErrY()==0) continue;
5195 //if (point->GetErrZ()==0) continue;
5197 //Float_t ddy = (point->GetY()-cl->GetY())/point->GetErrY();
5198 //Float_t ddz = (point->GetZ()-cl->GetZ())/point->GetErrZ();
5199 //if ((ddy*ddy+ddz*ddz)>10) continue;
5202 // if (point->GetCPoint().GetMax()<5) continue;
5203 if (cl->GetMax()<5) continue;
5204 Float_t angley = point->GetAngleY();
5205 Float_t anglez = point->GetAngleZ();
5207 Float_t rsigmay2 = point->GetSigmaY();
5208 Float_t rsigmaz2 = point->GetSigmaZ();
5212 rsigmay += pointm->GetTPoint().GetSigmaY();
5213 rsigmaz += pointm->GetTPoint().GetSigmaZ();
5217 rsigmay += pointp->GetTPoint().GetSigmaY();
5218 rsigmaz += pointp->GetTPoint().GetSigmaZ();
5225 Float_t rsigma = TMath::Sqrt(rsigmay2*rsigmaz2);
5227 Float_t ampc = 0; // normalization to the number of electrons
5229 // ampc = 1.*point->GetCPoint().GetMax();
5230 ampc = 1.*cl->GetMax();
5231 //ampc = 1.*point->GetCPoint().GetQ();
5232 // AliTPCClusterPoint & p = point->GetCPoint();
5233 // Float_t dy = TMath::Abs(Int_t( TMath::Abs(p.GetY()/0.6)) - TMath::Abs(p.GetY()/0.6)+0.5);
5234 // Float_t iz = (250.0-TMath::Abs(p.GetZ())+0.11)/0.566;
5236 // TMath::Abs( Int_t(iz) - iz + 0.5);
5237 //ampc *= 1.15*(1-0.3*dy);
5238 //ampc *= 1.15*(1-0.3*dz);
5239 // Float_t zfactor = (1.05-0.0004*TMath::Abs(point->GetCPoint().GetZ()));
5243 //ampc = 1.0*point->GetCPoint().GetMax();
5244 ampc = 1.0*cl->GetMax();
5245 //ampc = 1.0*point->GetCPoint().GetQ();
5246 //AliTPCClusterPoint & p = point->GetCPoint();
5247 // Float_t dy = TMath::Abs(Int_t( TMath::Abs(p.GetY()/0.4)) - TMath::Abs(p.GetY()/0.4)+0.5);
5248 //Float_t iz = (250.0-TMath::Abs(p.GetZ())+0.11)/0.566;
5250 // TMath::Abs( Int_t(iz) - iz + 0.5);
5252 //ampc *= 1.15*(1-0.3*dy);
5253 //ampc *= 1.15*(1-0.3*dz);
5254 // Float_t zfactor = (1.02-0.000*TMath::Abs(point->GetCPoint().GetZ()));
5258 ampc *= 2.0; // put mean value to channel 50
5259 //ampc *= 0.58; // put mean value to channel 50
5261 // if (type>0) w = 1./(type/2.-0.5);
5262 // Float_t z = TMath::Abs(cl->GetZ());
5265 //ampc /= (1+0.0008*z);
5269 //ampc /= (1+0.0008*z);
5271 //ampc /= (1+0.0008*z);
5274 if (type<0) { //amp at the border - lower weight
5279 if (rsigma>1.5) ampc/=1.3; // if big backround
5281 angular[nc[of]] = TMath::Sqrt(1.+angley*angley+anglez*anglez);
5286 TMath::Sort(nc[of],amp,index,kFALSE);
5290 //meanlog = amp[index[Int_t(nc[of]*0.33)]];
5292 for (Int_t i=int(nc[of]*low+0.5);i<int(nc[of]*up+0.5);i++){
5293 Float_t ampl = amp[index[i]]/angular[index[i]];
5294 ampl = meanlog*TMath::Log(1.+ampl/meanlog);
5296 sumw += weight[index[i]];
5297 sumamp += weight[index[i]]*ampl;
5298 sumamp2 += weight[index[i]]*ampl*ampl;
5299 norm[of] += angular[index[i]]*weight[index[i]];
5306 mean[of] = sumamp/sumw;
5307 sigma[of] = sumamp2/sumw-mean[of]*mean[of];
5309 sigma[of] = TMath::Sqrt(sigma[of]);
5313 mean[of] = (TMath::Exp(mean[of]/meanlog)-1)*meanlog;
5314 //mean *=(1-0.02*(sigma/(mean*0.17)-1.));
5315 //mean *=(1-0.1*(norm-1.));
5322 // mean[0]*= (1-0.05*(sigma[0]/(0.01+mean[1]*0.18)-1));
5323 // mean[1]*= (1-0.05*(sigma[1]/(0.01+mean[0]*0.18)-1));
5326 // dedx = (mean[0]* TMath::Sqrt((1.+nc[0]))+ mean[1]* TMath::Sqrt((1.+nc[1])) )/
5327 // ( TMath::Sqrt((1.+nc[0]))+TMath::Sqrt((1.+nc[1])));
5331 for (Int_t i =0;i<4;i++){
5332 if (nc[i]>2&&nc[i]<1000){
5333 dedx += mean[i] *nc[i];
5334 fSdEdx += sigma[i]*(nc[i]-2);
5335 fMAngular += norm[i] *nc[i];
5340 fSDEDX[i] = sigma[i];
5353 // Float_t dedx1 =dedx;
5356 for (Int_t i =0;i<4;i++){
5357 if (nc[i]>2&&nc[i]<1000){
5358 mean[i] = mean[i]*(1-0.12*(sigma[i]/(fSdEdx)-1.));
5359 dedx += mean[i] *nc[i];
5374 Double_t p=TMath::Sqrt((1.+ GetTgl()*GetTgl())/(Get1Pt()*Get1Pt()));
5377 if (dedx < 39.+ 12./(p+0.25)/(p+0.25)) { SetMass(0.13957); return;}
5378 if (dedx < 39.+ 12./p/p) { SetMass(0.49368); return;}
5379 SetMass(0.93827); return;
5383 if (dedx < 39.+ 12./(p+0.25)/(p+0.25)) { SetMass(0.13957); return;}
5384 SetMass(0.93827); return;
5387 SetMass(0.13957); return;
5397 void AliTPCseed::CookdEdx2(Double_t low, Double_t up) {
5398 //-----------------------------------------------------------------
5399 // This funtion calculates dE/dX within the "low" and "up" cuts.
5400 //-----------------------------------------------------------------
5403 Float_t angular[200];
5404 Float_t weight[200];
5406 Bool_t inlimit[200];
5407 for (Int_t i=0;i<200;i++) inlimit[i]=kFALSE;
5408 for (Int_t i=0;i<200;i++) amp[i]=10000;
5409 for (Int_t i=0;i<200;i++) angular[i]= 1;;
5413 Float_t meanlog = 100.;
5414 Int_t indexde[4]={0,64,128,160};
5421 Float_t mean[4] = {0,0,0,0};
5422 Float_t sigma[4] = {1000,1000,1000,1000};
5423 Int_t nc[4] = {0,0,0,0};
5424 Float_t norm[4] = {1000,1000,1000,1000};
5429 // for (Int_t of =0; of<3; of++){
5430 // for (Int_t i=indexde[of];i<indexde[of+1];i++)
5431 for (Int_t i =0; i<160;i++)
5433 AliTPCTrackPoint * point = GetTrackPoint(i);
5434 if (point==0) continue;
5435 if (point->fIsShared){
5439 Int_t type = point->GetCPoint().GetType();
5440 if (type<0) continue;
5441 if (point->GetCPoint().GetMax()<5) continue;
5442 Float_t angley = point->GetTPoint().GetAngleY();
5443 Float_t anglez = point->GetTPoint().GetAngleZ();
5444 Float_t rsigmay = point->GetCPoint().GetSigmaY();
5445 Float_t rsigmaz = point->GetCPoint().GetSigmaZ();
5446 Float_t rsigma = TMath::Sqrt(rsigmay*rsigmaz);
5448 Float_t ampc = 0; // normalization to the number of electrons
5450 ampc = point->GetCPoint().GetMax();
5453 ampc = point->GetCPoint().GetMax();
5455 ampc *= 2.0; // put mean value to channel 50
5456 // ampc *= 0.565; // put mean value to channel 50
5459 Float_t z = TMath::Abs(point->GetCPoint().GetZ());
5466 if (type<0) { //amp at the border - lower weight
5469 if (rsigma>1.5) ampc/=1.3; // if big backround
5470 angular[i] = TMath::Sqrt(1.+angley*angley+anglez*anglez);
5471 amp[i] = ampc/angular[i];
5476 TMath::Sort(159,amp,index,kFALSE);
5477 for (Int_t i=int(anc*low+0.5);i<int(anc*up+0.5);i++){
5478 inlimit[index[i]] = kTRUE; // take all clusters
5481 // meanlog = amp[index[Int_t(anc*0.3)]];
5483 for (Int_t of =0; of<3; of++){
5487 for (Int_t i=indexde[of];i<indexde[of+1];i++)
5489 if (inlimit[i]==kFALSE) continue;
5490 Float_t ampl = amp[i];
5492 ampl = meanlog*TMath::Log(1.+ampl/meanlog);
5495 sumamp += weight[i]*ampl;
5496 sumamp2 += weight[i]*ampl*ampl;
5497 norm[of] += angular[i]*weight[i];
5505 mean[of] = sumamp/sumw;
5506 sigma[of] = sumamp2/sumw-mean[of]*mean[of];
5508 sigma[of] = TMath::Sqrt(sigma[of]);
5511 mean[of] = (TMath::Exp(mean[of]/meanlog)-1)*meanlog;
5521 Float_t www[3] = {12.,14.,17.};
5522 //Float_t www[3] = {1.,1.,1.};
5524 for (Int_t i =0;i<3;i++){
5525 if (nc[i]>2&&nc[i]<1000){
5526 dedx += mean[i] *nc[i]*www[i]/sigma[i];
5527 fSdEdx += sigma[i]*(nc[i]-2)*www[i]/sigma[i];
5528 fMAngular += norm[i] *nc[i];
5529 norm2 += nc[i]*www[i]/sigma[i];
5530 norm3 += (nc[i]-2)*www[i]/sigma[i];
5533 fSDEDX[i] = sigma[i];
5546 // Float_t dedx1 =dedx;
5550 for (Int_t i =0;i<3;i++){
5551 if (nc[i]>2&&nc[i]<1000&&sigma[i]>3){
5552 //mean[i] = mean[i]*(1+0.08*(sigma[i]/(fSdEdx)-1.));
5553 dedx += mean[i] *(nc[i])/(sigma[i]);
5554 norm4 += (nc[i])/(sigma[i]);
5558 if (norm4>0) dedx /= norm4;