]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPC.cxx
Accept only signals > 0
[u/mrichter/AliRoot.git] / TPC / AliTPC.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.15  1999/10/08 06:26:53  fca
19 Removed ClustersIndex - not used anymore
20
21 Revision 1.14  1999/09/29 09:24:33  fca
22 Introduction of the Copyright and cvs Log
23
24 */
25
26 ///////////////////////////////////////////////////////////////////////////////
27 //                                                                           //
28 //  Time Projection Chamber                                                  //
29 //  This class contains the basic functions for the Time Projection Chamber  //
30 //  detector. Functions specific to one particular geometry are              //
31 //  contained in the derived classes                                         //
32 //                                                                           //
33 //Begin_Html
34 /*
35 <img src="picts/AliTPCClass.gif">
36 */
37 //End_Html
38 //                                                                           //
39 //                                                                          //
40 ///////////////////////////////////////////////////////////////////////////////
41
42 #include <TMath.h>
43 #include <TRandom.h>
44 #include <TVector.h>
45 #include <TMatrix.h>
46 #include <TGeometry.h>
47 #include <TNode.h>
48 #include <TTUBS.h>
49 #include <TObjectTable.h>
50 #include "TParticle.h"
51 #include "AliTPC.h"
52 #include "AliRun.h"
53 #include <iostream.h>
54 #include <fstream.h>
55 #include "AliMC.h"
56
57 //MI change
58 #include "AliTPCParam.h"
59 #include "AliTPCD.h"
60 #include "AliTPCPRF2D.h"
61 #include "AliTPCRF1D.h"
62
63
64
65 ClassImp(AliTPC) 
66
67 //_____________________________________________________________________________
68 AliTPC::AliTPC()
69 {
70   //
71   // Default constructor
72   //
73   fIshunt   = 0;
74   fClusters = 0;
75   fHits     = 0;
76   fDigits   = 0;
77   fTracks   = 0;
78   fNsectors = 0;
79   fNtracks  = 0;
80   fNclusters= 0;
81
82   fDigParam= new AliTPCD();
83   fDigits = fDigParam->GetArray();
84 }
85  
86 //_____________________________________________________________________________
87 AliTPC::AliTPC(const char *name, const char *title)
88       : AliDetector(name,title)
89 {
90   //
91   // Standard constructor
92   //
93
94   //
95   // Initialise arrays of hits and digits 
96   fHits     = new TClonesArray("AliTPChit",  176);
97   //  fDigits   = new TClonesArray("AliTPCdigit",10000);
98   //MI change
99   fDigParam= new AliTPCD;
100   fDigits = fDigParam->GetArray();
101
102   AliTPCParam  *fTPCParam = &(fDigParam->GetParam());
103
104   //
105   // Initialise counters
106   //
107   fClusters = 0;
108   fTracks   = 0;
109   fNsectors = fTPCParam->GetNSector();
110   fNtracks  = 0;
111   fNclusters= 0;
112   fDigitsIndex = new Int_t[fNsectors+1];
113   //
114   fIshunt     =  0;
115   //
116   // Initialise color attributes
117   SetMarkerColor(kYellow);
118 }
119
120 //_____________________________________________________________________________
121 AliTPC::~AliTPC()
122 {
123   //
124   // TPC destructor
125   //
126   fIshunt   = 0;
127   delete fHits;
128   delete fDigits;
129   delete fClusters;
130   delete fTracks;
131   delete fDigParam;
132   if (fDigitsIndex)   delete [] fDigitsIndex;
133 }
134
135 //_____________________________________________________________________________
136 void AliTPC::AddCluster(Float_t *hits, Int_t *tracks)
137 {
138   //
139   // Add a simulated cluster to the list
140   //
141   if(!fClusters) fClusters=new TClonesArray("AliTPCcluster",10000);
142   TClonesArray &lclusters = *fClusters;
143   new(lclusters[fNclusters++]) AliTPCcluster(hits,tracks);
144 }
145  
146 //_____________________________________________________________________________
147 void AliTPC::AddCluster(const AliTPCcluster &c)
148 {
149   //
150   // Add a simulated cluster copy to the list
151   //
152   if(!fClusters) fClusters=new TClonesArray("AliTPCcluster",10000);
153   TClonesArray &lclusters = *fClusters;
154   new(lclusters[fNclusters++]) AliTPCcluster(c);
155 }
156  
157 //_____________________________________________________________________________
158 void AliTPC::AddDigit(Int_t *tracks, Int_t *digits)
159 {
160   //
161   // Add a TPC digit to the list
162   //
163   //  TClonesArray &ldigits = *fDigits;
164   //MI change 
165   TClonesArray &ldigits = *fDigParam->GetArray();
166   new(ldigits[fNdigits++]) AliTPCdigit(tracks,digits);
167 }
168  
169 //_____________________________________________________________________________
170 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
171 {
172   //
173   // Add a hit to the list
174   //
175   TClonesArray &lhits = *fHits;
176   new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
177 }
178  
179 //_____________________________________________________________________________
180 void AliTPC::AddTrack(Float_t *hits)
181 {
182   //
183   // Add a track to the list of tracks
184   //
185   TClonesArray &ltracks = *fTracks;
186   new(ltracks[fNtracks++]) AliTPCtrack(hits);
187 }
188
189 //_____________________________________________________________________________
190 void AliTPC::AddTrack(const AliTPCtrack& t)
191 {
192   //
193   // Add a track copy to the list of tracks
194   //
195   if(!fTracks) fTracks=new TClonesArray("AliTPCtrack",10000);
196   TClonesArray &ltracks = *fTracks;
197   new(ltracks[fNtracks++]) AliTPCtrack(t);
198 }
199
200 //_____________________________________________________________________________
201 void AliTPC::BuildGeometry()
202 {
203   //
204   // Build TPC ROOT TNode geometry for the event display
205   //
206   TNode *Node, *Top;
207   TTUBS *tubs;
208   Int_t i;
209   const int kColorTPC=19;
210   char name[5], title[25];
211   const Double_t kDegrad=TMath::Pi()/180;
212   const Double_t kRaddeg=180./TMath::Pi();
213
214   AliTPCParam * fTPCParam = &(fDigParam->GetParam());
215
216   Float_t InnerOpenAngle = fTPCParam->GetInnerAngle();
217   Float_t OuterOpenAngle = fTPCParam->GetOuterAngle();
218
219   Float_t InnerAngleShift = fTPCParam->GetInnerAngleShift();
220   Float_t OuterAngleShift = fTPCParam->GetOuterAngleShift();
221
222   Int_t nLo = fTPCParam->GetNInnerSector()/2;
223   Int_t nHi = fTPCParam->GetNOuterSector()/2;  
224
225   const Double_t loAng = (Double_t)TMath::Nint(InnerOpenAngle*kRaddeg);
226   const Double_t hiAng = (Double_t)TMath::Nint(OuterOpenAngle*kRaddeg);
227   const Double_t loAngSh = (Double_t)TMath::Nint(InnerAngleShift*kRaddeg);
228   const Double_t hiAngSh = (Double_t)TMath::Nint(OuterAngleShift*kRaddeg);  
229
230
231   const Double_t loCorr = 1/TMath::Cos(0.5*loAng*kDegrad);
232   const Double_t hiCorr = 1/TMath::Cos(0.5*hiAng*kDegrad);
233
234   Double_t rl,ru;
235   
236
237   //
238   // Get ALICE top node
239   //
240
241   Top=gAlice->GetGeometry()->GetNode("alice");
242
243   //  inner sectors
244
245   rl = fTPCParam->GetInSecLowEdge();
246   ru = fTPCParam->GetInSecUpEdge();
247  
248
249   for(i=0;i<nLo;i++) {
250     sprintf(name,"LS%2.2d",i);
251     name[4]='\0';
252     sprintf(title,"TPC low sector %3d",i);
253     title[24]='\0';
254     
255     tubs = new TTUBS(name,title,"void",rl*loCorr,ru*loCorr,250.,
256                      loAng*(i-0.5)+loAngSh,loAng*(i+0.5)+loAngSh);
257     tubs->SetNumberOfDivisions(1);
258     Top->cd();
259     Node = new TNode(name,title,name,0,0,0,"");
260     Node->SetLineColor(kColorTPC);
261     fNodes->Add(Node);
262   }
263
264   // Outer sectors
265
266   rl = fTPCParam->GetOuSecLowEdge();
267   ru = fTPCParam->GetOuSecUpEdge();
268
269   for(i=0;i<nHi;i++) {
270     sprintf(name,"US%2.2d",i);
271     name[4]='\0';
272     sprintf(title,"TPC upper sector %d",i);
273     title[24]='\0';
274     tubs = new TTUBS(name,title,"void",rl*hiCorr,ru*hiCorr,250,
275                      hiAng*(i-0.5)+hiAngSh,hiAng*(i+0.5)+hiAngSh);
276     tubs->SetNumberOfDivisions(1);
277     Top->cd();
278     Node = new TNode(name,title,name,0,0,0,"");
279     Node->SetLineColor(kColorTPC);
280     fNodes->Add(Node);
281   }
282 }  
283   
284   
285
286 //_____________________________________________________________________________
287 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
288 {
289   //
290   // Calculate distance from TPC to mouse on the display
291   // Dummy procedure
292   //
293   return 9999;
294 }
295
296 //_____________________________________________________________________________
297 static Double_t SigmaY2(Double_t r, Double_t tgl, Double_t pt)
298 {
299   //
300   // Parametrised error of the cluster reconstruction (pad direction)   
301   //
302   pt=TMath::Abs(pt)*1000.;
303   Double_t x=r/pt;
304   tgl=TMath::Abs(tgl);
305   Double_t s=a_rphi - b_rphi*r*tgl + c_rphi*x*x + d_rphi*x;
306   if (s<0.4e-3) s=0.4e-3;
307   s*=1.3; //Iouri Belikov
308   return s;
309 }
310
311 //_____________________________________________________________________________
312 static Double_t SigmaZ2(Double_t r, Double_t tgl) 
313 {
314   //
315   // Parametrised error of the cluster reconstruction (drift direction)
316   //
317   tgl=TMath::Abs(tgl);
318   Double_t s=a_z - b_z*r*tgl + c_z*tgl*tgl;
319   if (s<0.4e-3) s=0.4e-3;
320   s*=1.3; //Iouri Belikov
321   return s;
322 }
323
324 //_____________________________________________________________________________
325 inline Double_t f1(Double_t x1,Double_t y1,
326                    Double_t x2,Double_t y2,
327                    Double_t x3,Double_t y3) 
328 {
329   //-----------------------------------------------------------------
330   // Initial approximation of the track curvature
331   //
332   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
333   //-----------------------------------------------------------------
334   Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
335   Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
336                   (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
337   Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
338                   (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
339
340   Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
341
342   return -xr*yr/sqrt(xr*xr+yr*yr); 
343 }
344
345
346 //_____________________________________________________________________________
347 inline Double_t f2(Double_t x1,Double_t y1,
348                    Double_t x2,Double_t y2,
349                    Double_t x3,Double_t y3) 
350 {
351   //-----------------------------------------------------------------
352   // Initial approximation of the track curvature times center of curvature
353   //
354   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
355   //-----------------------------------------------------------------
356   Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
357   Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
358                   (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
359   Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
360                   (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
361
362   Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
363   
364   return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
365 }
366
367 //_____________________________________________________________________________
368 inline Double_t f3(Double_t x1,Double_t y1, 
369                    Double_t x2,Double_t y2,
370                    Double_t z1,Double_t z2) 
371 {
372   //-----------------------------------------------------------------
373   // Initial approximation of the tangent of the track dip angle
374   //
375   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
376   //-----------------------------------------------------------------
377   return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
378 }
379
380 //_____________________________________________________________________________
381 static int FindProlongation(AliTPCtrack& t, const AliTPCSector *sec,
382                             int s, int rf=0) 
383 {
384   //-----------------------------------------------------------------
385   // This function tries to find a track prolongation.
386   //
387   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
388   //-----------------------------------------------------------------
389   const int ROWS_TO_SKIP=int(0.5*sec->GetNRows());
390   const Float_t MAX_CHI2=12.;
391   int try_again=ROWS_TO_SKIP;
392   Double_t alpha=sec->GetAlpha();
393   int ns=int(2*TMath::Pi()/alpha+0.5);
394
395   for (int nr=sec->GetRowNumber(t.GetX())-1; nr>=rf; nr--) {
396     Double_t x=sec->GetX(nr), ymax=sec->GetMaxY(nr);
397     if (!t.PropagateTo(x)) return 0;
398
399     AliTPCcluster *cl=0;
400     Double_t max_chi2=MAX_CHI2;
401     const AliTPCRow& row=sec[s][nr];
402     Double_t sy2=SigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
403     Double_t sz2=SigmaZ2(t.GetX(),t.GetTgl());
404     Double_t road=5.*sqrt(t.GetSigmaY2() + sy2), y=t.GetY(), z=t.GetZ();
405
406     if (road>30) {
407       if (t>4) cerr<<t<<" FindProlongation warning: Too broad road !\n"; 
408       return 0;
409     }
410
411     if (row) {
412       for (int i=row.Find(y-road); i<row; i++) {
413         AliTPCcluster* c=(AliTPCcluster*)(row[i]);
414         if (c->fY > y+road) break;
415         if (c->IsUsed()) continue;
416         if ((c->fZ - z)*(c->fZ - z) > 25.*(t.GetSigmaZ2() + sz2)) continue;
417         Double_t chi2=t.GetPredictedChi2(c);
418         if (chi2 > max_chi2) continue;
419         max_chi2=chi2;
420         cl=c;       
421       }
422     }
423     if (cl) {
424       t.Update(cl,max_chi2);
425       Double_t ll=TMath::Sqrt((1+t.GetTgl()*t.GetTgl())/
426                (1-(t.GetC()*x-t.GetEta())*(t.GetC()*x-t.GetEta())));
427       cl->fdEdX = cl->fQ/ll;
428       try_again=ROWS_TO_SKIP;
429     } else {
430       if (try_again==0) break;
431       if (y > ymax) {
432          s = (s+1) % ns;
433          if (!t.Rotate(alpha)) return 0;
434       } else if (y <-ymax) {
435          s = (s-1+ns) % ns;
436          if (!t.Rotate(-alpha)) return 0;
437       }
438       try_again--;
439     }
440   }
441
442   return 1;
443 }
444
445
446 //_____________________________________________________________________________
447 static void MakeSeeds(TObjArray& seeds,const AliTPCSector *sec, int max_sec,
448 int i1, int i2)
449 {
450   //-----------------------------------------------------------------
451   // This function creates track seeds.
452   //
453   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
454   //-----------------------------------------------------------------
455   TMatrix C(5,5); TVector x(5);
456   double alpha=sec->GetAlpha(), shift=sec->GetAlphaShift();
457   double cs=cos(alpha), sn=sin(alpha);
458   for (int ns=0; ns<max_sec; ns++) {
459     int nl=sec[(ns-1+max_sec)%max_sec][i2];
460     int nm=sec[ns][i2];
461     int nu=sec[(ns+1)%max_sec][i2];
462     const AliTPCRow& r1=sec[ns][i1];
463     for (int is=0; is < r1; is++) {
464       double x1=sec->GetX(i1), y1=r1[is]->fY, z1=r1[is]->fZ;
465       for (int js=0; js < nl+nm+nu; js++) {
466         const AliTPCcluster *cl;
467         int ks;
468         double x2=sec->GetX(i2), y2, z2, tmp;
469
470         if (js<nl) {
471           ks=(ns-1+max_sec)%max_sec;
472           const AliTPCRow& r2=sec[(ns-1+max_sec)%max_sec][i2];
473           cl=r2[js];
474           y2=cl->fY; z2=cl->fZ;
475           tmp= x2*cs+y2*sn;
476           y2 =-x2*sn+y2*cs; x2=tmp;
477         } else 
478           if (js<nl+nm) {
479             ks=ns;
480             const AliTPCRow& r2=sec[ns][i2];
481             cl=r2[js-nl];
482             y2=cl->fY; z2=cl->fZ;
483           } else {
484             ks=(ns+1)%max_sec;
485             const AliTPCRow& r2=sec[(ns+1)%max_sec][i2];
486             cl=r2[js-nl-nm];
487             y2=cl->fY; z2=cl->fZ;
488             tmp=x2*cs-y2*sn;
489             y2 =x2*sn+y2*cs; x2=tmp;
490           }
491
492         double d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
493         if (d==0.) {cerr<<"MakeSeeds warning: Straight seed !\n"; continue;}
494
495         x(0)=y1;
496         x(1)=z1;
497         x(2)=f1(x1,y1,x2,y2,0.,0.);
498         x(3)=f2(x1,y1,x2,y2,0.,0.);
499         x(4)=f3(x1,y1,x2,y2,z1,z2);
500         
501         if (TMath::Abs(x(2)*x1-x(3)) >= 0.999) continue;
502         
503         if (TMath::Abs(x(4)) > 1.2) continue;
504
505         Double_t a=asin(x(3));
506         Double_t zv=z1 - x(4)/x(2)*(a+asin(x(2)*x1-x(3)));
507         if (TMath::Abs(zv)>33.) continue; 
508
509         TMatrix X(6,6); X=0.; 
510         X(0,0)=r1[is]->fSigmaY2; X(1,1)=r1[is]->fSigmaZ2;
511         X(2,2)=cl->fSigmaY2;     X(3,3)=cl->fSigmaZ2;
512         X(4,4)=3./12.; X(5,5)=3./12.;
513         TMatrix F(5,6); F.UnitMatrix();
514         Double_t sy=sqrt(X(0,0)), sz=sqrt(X(1,1));
515         F(2,0)=(f1(x1,y1+sy,x2,y2,0.,0.)-x(2))/sy;
516         F(2,2)=(f1(x1,y1,x2,y2+sy,0.,0.)-x(2))/sy;
517         F(2,4)=(f1(x1,y1,x2,y2,0.,0.+sy)-x(2))/sy;
518         F(3,0)=(f2(x1,y1+sy,x2,y2,0.,0.)-x(3))/sy;
519         F(3,2)=(f2(x1,y1,x2,y2+sy,0.,0.)-x(3))/sy;
520         F(3,4)=(f2(x1,y1,x2,y2,0.,0.+sy)-x(3))/sy;
521         F(4,0)=(f3(x1,y1+sy,x2,y2,z1,z2)-x(4))/sy;
522         F(4,1)=(f3(x1,y1,x2,y2,z1+sz,z2)-x(4))/sz;
523         F(4,2)=(f3(x1,y1,x2,y2+sy,z1,z2)-x(4))/sy;
524         F(4,3)=(f3(x1,y1,x2,y2,z1,z2+sz)-x(4))/sz;
525         F(4,4)=0;
526         F(3,3)=0;
527         
528         TMatrix t(F,TMatrix::kMult,X);
529         C.Mult(t,TMatrix(TMatrix::kTransposed,F));
530
531         AliTPCtrack *track=new AliTPCtrack(r1[is], x, C, x1, ns*alpha+shift);
532         int rc=FindProlongation(*track,sec,ns,i2);
533         if (rc<0 || *track<(i1-i2)/2) delete track;
534         else seeds.AddLast(track); 
535       }
536     }
537   }
538 }
539
540 //_____________________________________________________________________________
541 AliTPCParam *AliTPCSector::param;
542 void AliTPC::Clusters2Tracks()
543 {
544   //-----------------------------------------------------------------
545   // This is a track finder.
546   //
547   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
548   //-----------------------------------------------------------------
549   if (!fClusters) return;
550
551   AliTPCParam *p=&fDigParam->GetParam();
552   AliTPCSector::SetParam(p);
553
554   const int nis=p->GetNInnerSector()/2;
555   AliTPCSSector *ssec=new AliTPCSSector[nis];         
556   int nrow_low=ssec->GetNRows();     
557
558   const int nos=p->GetNOuterSector()/2;
559   AliTPCLSector *lsec=new AliTPCLSector[nos];
560   int nrow_up=lsec->GetNRows();
561
562   int ncl=fClusters->GetEntriesFast();
563   while (ncl--) {
564     AliTPCcluster *c=(AliTPCcluster*)fClusters->UncheckedAt(ncl);
565     Int_t sec=c->fSector, row=c->fPadRow;
566     if (sec<nis*2) {
567       ssec[sec%nis][row].InsertCluster(c);
568     } else {
569       sec -= nis*2;
570       lsec[sec%nos][row].InsertCluster(c);
571     }
572   }
573   
574   TObjArray seeds(20000);
575
576   int nrows=nrow_low+nrow_up;
577   int gap=int(0.125*nrows), shift=int(0.5*gap);
578   MakeSeeds(seeds, lsec, nos, nrow_up-1, nrow_up-1-gap);
579   MakeSeeds(seeds, lsec, nos, nrow_up-1-shift, nrow_up-1-shift-gap);
580     
581   seeds.Sort();
582   
583   int found=0;
584   int nseed=seeds.GetEntriesFast();
585   
586   for (int s=0; s<nseed; s++) {
587     AliTPCtrack& t=*((AliTPCtrack*)seeds.UncheckedAt(s));
588     double alpha=t.GetAlpha();
589     if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();  
590     if (alpha < 0.            ) alpha += 2.*TMath::Pi();  
591     int ns=int(alpha/lsec->GetAlpha())%nos;
592
593     if (!FindProlongation(t,lsec,ns)) continue;
594
595     alpha=t.GetAlpha() + 0.5*ssec->GetAlpha() - ssec->GetAlphaShift();
596     if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
597     if (alpha < 0.            ) alpha += 2.*TMath::Pi();
598     ns=int(alpha/ssec->GetAlpha())%nis; //index of the inner sector needed
599
600     alpha=ns*ssec->GetAlpha() - t.GetAlpha();
601     if (!t.Rotate(alpha)) continue;
602
603     if (!FindProlongation(t,ssec,ns)) continue;
604     
605     if (t < int(0.4*nrows)) continue;
606     
607     AddTrack(t);
608     t.UseClusters();
609     cerr<<found++<<'\r';
610   }  
611
612   delete[] ssec;
613   delete[] lsec;
614 }
615
616 //_____________________________________________________________________________
617 void AliTPC::CreateMaterials()
618 {
619   //-----------------------------------------------
620   // Create Materials for for TPC
621   //-----------------------------------------------
622
623   //-----------------------------------------------------------------
624   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
625   //-----------------------------------------------------------------
626
627   Int_t ISXFLD=gAlice->Field()->Integ();
628   Float_t SXMGMX=gAlice->Field()->Max();
629
630   Float_t amat[5]; // atomic numbers
631   Float_t zmat[5]; // z
632   Float_t wmat[5]; // proportions
633
634   Float_t density;
635
636   //  ********************* Gases *******************
637
638   //--------------------------------------------------------------
639   // pure gases
640   //--------------------------------------------------------------
641
642   // Ne
643
644
645   Float_t a_ne = 20.18;
646   Float_t z_ne = 10.;
647   
648   density = 0.0009;
649
650   AliMaterial(20,"Ne",a_ne,z_ne,density,999.,999.);
651
652   // Ar
653
654   Float_t a_ar = 39.948;
655   Float_t z_ar = 18.;
656
657   density = 0.001782;
658  
659   AliMaterial(21,"Ar",a_ar,z_ar,density,999.,999.);
660
661   Float_t a_pure[2];
662   
663   a_pure[0] = a_ne;
664   a_pure[1] = a_ar;
665   
666
667   //--------------------------------------------------------------
668   // gases - compounds
669   //--------------------------------------------------------------
670
671   Float_t amol[3];
672
673   //  CO2
674
675   amat[0]=12.011;
676   amat[1]=15.9994;
677
678   zmat[0]=6.;
679   zmat[1]=8.;
680
681   wmat[0]=1.;
682   wmat[1]=2.;
683
684   density=0.001977;
685
686   amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
687
688   AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
689
690   // CF4
691
692   amat[0]=12.011;
693   amat[1]=18.998;
694
695   zmat[0]=6.;
696   zmat[1]=9.;
697  
698   wmat[0]=1.;
699   wmat[1]=4.;
700  
701   density=0.003034;
702
703   amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
704
705   AliMixture(11,"CF4",amat,zmat,density,-2,wmat); 
706
707   // CH4
708
709   amat[0]=12.011;
710   amat[1]=1.;
711
712   zmat[0]=6.;
713   zmat[1]=1.;
714
715   wmat[0]=1.;
716   wmat[1]=4.;
717
718   density=0.000717;
719
720   amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
721
722   AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
723
724   //----------------------------------------------------------------
725   // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
726   //----------------------------------------------------------------
727  
728   char namate[21];
729  
730   density = 0.;
731   Float_t am=0;
732   Int_t nc;
733
734   Float_t a,z,rho,absl,X0,buf[1];
735   Int_t nbuf;
736
737   for(nc = 0;nc<fNoComp;nc++)
738     {
739     
740       // retrive material constants
741       
742       gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
743
744       amat[nc] = a;
745       zmat[nc] = z;
746
747       Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
748  
749       am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? a_pure[nnc] : amol[nnc]); 
750       density += fMixtProp[nc]*rho;  // density of the mixture
751       
752     }
753
754   // mixture proportions by weight!
755
756   for(nc = 0;nc<fNoComp;nc++)
757     {
758
759       Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
760
761       wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ? a_pure[nnc] : amol[nnc])/am;
762
763     }  
764   
765   AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
766   AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
767   AliMixture(33,"Drift gas 3",amat,zmat,density,fNoComp,wmat); 
768
769   AliMedium(2, "Drift gas 1", 31, 0, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
770   AliMedium(3, "Drift gas 2", 32, 0, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
771   AliMedium(4, "Drift gas 3", 33, 1, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
772
773   // Air 
774
775   AliMaterial(24, "Air", 14.61, 7.3, .001205, 30420., 67500.);
776
777   AliMedium(24, "Air", 24, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
778
779   //----------------------------------------------------------------------
780   //               solid materials
781   //----------------------------------------------------------------------
782
783   // Al
784
785   AliMaterial(30, "Al", 26.98, 13., 2.7, 8.9, 37.2);
786
787   AliMedium(0, "Al",30, 0, ISXFLD, SXMGMX, 10., .1, .1, .1,   .1);
788
789   // Si
790
791   AliMaterial(31, "Si", 28.086, 14.,2.33, 9.36, 999.);
792
793   AliMedium(7, "Al",31, 0, ISXFLD, SXMGMX, 10., .1, .1, .1,   .1);
794   
795
796   // Mylar C5H4O2
797
798   amat[0]=12.011;
799   amat[1]=1.;
800   amat[2]=15.9994;
801
802   zmat[0]=6.;
803   zmat[1]=1.;
804   zmat[2]=8.;
805
806   wmat[0]=5.;
807   wmat[1]=4.;
808   wmat[2]=2.; 
809
810   density = 1.39;
811   
812   AliMixture(32, "Mylar",amat,zmat,density,-3,wmat);
813
814   AliMedium(5, "Mylar",32, 0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
815
816
817
818
819   // Carbon (normal)
820
821   AliMaterial(33,"C normal",12.011,6.,2.265,18.8,999.);
822
823   AliMedium(6,"C normal",33,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
824
825   // G10 for inner and outr field cage
826   // G10 is 60% SiO2 + 40% epoxy, right now I use A and Z for SiO2
827
828   Float_t rhoFactor;
829
830   amat[0]=28.086;
831   amat[1]=15.9994;
832
833   zmat[0]=14.;
834   zmat[1]=8.;
835
836   wmat[0]=1.;
837   wmat[1]=2.;
838
839   density = 1.7;
840   
841
842   AliMixture(34,"G10 aux.",amat,zmat,density,-2,wmat);
843
844
845   gMC->Gfmate((*fIdmate)[34],namate,a,z,rho,X0,absl,buf,nbuf);
846
847   Float_t thickX0 = 0.0052; // field cage in X0 units
848   
849   Float_t thick = 2.; // in cm
850
851   X0=19.4; // G10 
852
853   rhoFactor = X0*thickX0/thick;
854   density = rho*rhoFactor;
855
856   AliMaterial(35,"G10-fc",a,z,density,999.,999.);
857
858   AliMedium(8,"G10-fc",35,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
859
860   thickX0 = 0.0027; // inner vessel (eta <0.9)
861   thick=0.5;
862   rhoFactor = X0*thickX0/thick;
863   density = rho*rhoFactor;
864
865   AliMaterial(36,"G10-iv",a,z,density,999.,999.);  
866
867   AliMedium(9,"G10-iv",36,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
868
869   //  Carbon fibre  
870   
871   gMC->Gfmate((*fIdmate)[33],namate,a,z,rho,X0,absl,buf,nbuf);
872
873   thickX0 = 0.0133; // outer vessel
874   thick=3.0;
875   rhoFactor = X0*thickX0/thick;
876   density = rho*rhoFactor;
877
878
879   AliMaterial(37,"C-ov",a,z,density,999.,999.);
880
881   AliMedium(10,"C-ov",37,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);  
882
883   thickX0=0.015; // inner vessel (cone, eta > 0.9)
884   thick=1.5;
885   rhoFactor = X0*thickX0/thick;
886   density = rho*rhoFactor;
887
888   AliMaterial(38,"C-ivc",a,z,density,999.,999.);
889
890   AliMedium(11,"C-ivc",38,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
891
892   //
893
894   AliMedium(12,"CO2",10,0, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
895     
896
897
898 }
899
900 //_____________________________________________________________________________
901 struct Bin {
902    const AliTPCdigit *dig;
903    int idx;
904    Bin() {dig=0; idx=-1;}
905 };
906
907 struct PreCluster : public AliTPCcluster {
908   const AliTPCdigit* summit; //pointer to the maximum digit of this precluster
909   int idx;                   //index in AliTPC::fClusters
910   int npeaks;                //number of peaks in this precluster
911   int ndigits;               //number of digits in this precluster
912   PreCluster();
913 };
914 PreCluster::PreCluster() : AliTPCcluster() {npeaks=ndigits=0;}
915
916 //_____________________________________________________________________________
917 static void FindPreCluster(int i,int j,int maxj,Bin *bins,PreCluster &c) 
918 {
919   //-----------------------------------------------------------------
920   // This function looks for "preclusters".
921   //
922   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
923   //-----------------------------------------------------------------
924   Bin& b=bins[i*maxj+j];
925   double q=(double)TMath::Abs(b.dig->fSignal);
926
927   if (b.idx >= 0 && b.idx != c.idx) {
928     c.idx=b.idx;
929     c.npeaks++;
930   }
931   
932   if (q > TMath::Abs(c.summit->fSignal)) c.summit=b.dig;
933   
934   c.fY += i*q;
935   c.fZ += j*q;
936   c.fSigmaY2 += i*i*q;
937   c.fSigmaZ2 += j*j*q;
938   c.fQ += q;
939   c.ndigits++;
940
941   b.dig = 0;  b.idx = c.idx;
942   
943   if (bins[(i-1)*maxj+j].dig) FindPreCluster(i-1,j,maxj,bins,c);
944   if (bins[i*maxj+(j-1)].dig) FindPreCluster(i,j-1,maxj,bins,c);
945   if (bins[(i+1)*maxj+j].dig) FindPreCluster(i+1,j,maxj,bins,c);
946   if (bins[i*maxj+(j+1)].dig) FindPreCluster(i,j+1,maxj,bins,c);
947   
948 }
949
950 //_____________________________________________________________________________
951 void AliTPC::Digits2Clusters()
952 {
953   //-----------------------------------------------------------------
954   // This is a simple cluster finder.
955   //
956   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
957   //-----------------------------------------------------------------
958   AliTPCParam *par = &(fDigParam->GetParam());
959   
960   int inp=par->GetNPads(0,                  par->GetNRowLow()-1);
961   int onp=par->GetNPads(par->GetNSector()-1,par->GetNRowUp() -1);
962   const int MAXY=(inp>onp) ? inp+2 : onp+2;
963   const int MAXTBKT=int((z_end+6*par->GetZSigma())/par->GetZWidth())+1;
964   const int MAXZ=MAXTBKT+2;
965   const int THRESHOLD=20;
966   
967   TTree *t=(TTree*)gDirectory->Get("TreeD0_Param1");
968   t->GetBranch("Digits")->SetAddress(&fDigits);
969   Int_t sectors_by_rows=(Int_t)t->GetEntries();
970   
971   int ncls=0;
972   
973   for (Int_t n=0; n<sectors_by_rows; n++) {
974     if (!t->GetEvent(n)) continue;
975     Bin *bins=new Bin[MAXY*MAXZ];
976     AliTPCdigit *dig=(AliTPCdigit*)fDigits->UncheckedAt(0);
977     int sec=dig->fSector, row=dig->fPadRow;
978     int ndigits=fDigits->GetEntriesFast();
979     
980     int npads, sign;
981     {
982        int nis=par->GetNInnerSector(), nos=par->GetNOuterSector();
983        if (sec < nis) {
984           npads = par->GetNPadsLow(row);
985           sign = (sec < nis/2) ? 1 : -1;
986        } else {
987           npads = par->GetNPadsUp(row);
988           sign = ((sec-nis) < nos/2) ? 1 : -1;
989        }
990     }
991
992     int ndig;
993     for (ndig=0; ndig<ndigits; ndig++) {
994       dig=(AliTPCdigit*)fDigits->UncheckedAt(ndig);
995       int i=dig->fPad+1, j=dig->fTime+1;
996       if (i > npads) {
997          cerr<<"AliTPC::Digits2Clusters error: pad number is out of range ! ";
998          cerr<<i<<' '<<npads<<endl; 
999          continue;
1000       }
1001       if (j > MAXTBKT) {
1002          cerr<<"AliTPC::Digits2Clusters error: time bucket is out of range ! ";
1003          cerr<<j<<' '<<MAXTBKT<<endl; 
1004          continue;
1005       }
1006       if (dig->fSignal >= THRESHOLD) bins[i*MAXZ+j].dig=dig;
1007       if (i==1 || i==npads || j==1 || j==MAXTBKT) dig->fSignal*=-1;
1008     }
1009
1010     int ncl=0;
1011     int i,j;
1012     
1013     for (i=1; i<MAXY-1; i++) {
1014       for (j=1; j<MAXZ-1; j++) {
1015         if (bins[i*MAXZ+j].dig == 0) continue;
1016         PreCluster c; c.summit=bins[i*MAXZ+j].dig; c.idx=ncls;
1017         FindPreCluster(i, j, MAXZ, bins, c);
1018         c.fY /= c.fQ;
1019         c.fZ /= c.fQ;
1020
1021         double s2 = c.fSigmaY2/c.fQ - c.fY*c.fY;
1022         c.fSigmaY2 = s2 + 1./12.;
1023         c.fSigmaY2 *= par->GetPadPitchWidth()*par->GetPadPitchWidth();
1024         if (s2 != 0.) c.fSigmaY2 *= 0.17;
1025
1026         s2 = c.fSigmaZ2/c.fQ - c.fZ*c.fZ;
1027         c.fSigmaZ2 = s2 + 1./12.;
1028         c.fSigmaZ2 *= par->GetZWidth()*par->GetZWidth();
1029         if (s2 != 0.) c.fSigmaZ2 *= 0.41;
1030
1031         c.fY = (c.fY - 0.5 - 0.5*npads)*par->GetPadPitchWidth();
1032         c.fZ = par->GetZWidth()*c.fZ; 
1033         c.fZ -= 3.*par->GetZSigma(); // PASA delay 
1034         c.fZ = sign*(z_end - c.fZ);
1035
1036         c.fSector=sec;
1037         c.fPadRow=row;
1038         c.fTracks[0]=c.summit->fTracks[0];
1039         c.fTracks[1]=c.summit->fTracks[1];
1040         c.fTracks[2]=c.summit->fTracks[2];
1041
1042         if (c.summit->fSignal<0) {
1043           c.fSigmaY2 *= 25.;
1044           c.fSigmaZ2 *= 4.;
1045         }
1046         
1047         AddCluster(c); ncls++; ncl++;
1048       }
1049     }
1050         
1051     for (ndig=0; ndig<ndigits; ndig++) {
1052       dig=(AliTPCdigit*)fDigits->UncheckedAt(ndig);
1053       int i=dig->fPad+1, j=dig->fTime+1;
1054       if (i > npads) {
1055          cerr<<"AliTPC::Digits2Clusters error: pad number is out of range ! ";
1056          cerr<<i<<' '<<npads<<endl; 
1057          continue;
1058       }
1059       if (j > MAXTBKT) {
1060          cerr<<"AliTPC::Digits2Clusters error: time bucket is out of range ! ";
1061          cerr<<j<<' '<<MAXTBKT<<endl; 
1062          continue;
1063       }
1064       if (TMath::Abs(dig->fSignal)>=par->GetZeroSup()) bins[i*MAXZ+j].dig=dig;
1065     }
1066     
1067     for (i=1; i<MAXY-1; i++) {
1068       for (j=1; j<MAXZ-1; j++) {
1069         if (bins[i*MAXZ+j].dig == 0) continue;
1070         PreCluster c; c.summit=bins[i*MAXZ+j].dig; c.idx=ncls;
1071         FindPreCluster(i, j, MAXZ, bins, c);
1072         if (c.ndigits < 2) continue; //noise cluster
1073         if (c.npeaks>1) continue;    //overlapped cluster
1074         c.fY /= c.fQ;
1075         c.fZ /= c.fQ;
1076
1077         double s2 = c.fSigmaY2/c.fQ - c.fY*c.fY;
1078         c.fSigmaY2 = s2 + 1./12.;
1079         c.fSigmaY2 *= par->GetPadPitchWidth()*par->GetPadPitchWidth();
1080         if (s2 != 0.) c.fSigmaY2 *= 0.04;
1081
1082         s2 = c.fSigmaZ2/c.fQ - c.fZ*c.fZ;
1083         c.fSigmaZ2 = s2 + 1./12.;
1084         c.fSigmaZ2 *= par->GetZWidth()*par->GetZWidth();
1085         if (s2 != 0.) c.fSigmaZ2 *= 0.10;
1086
1087         c.fY = (c.fY - 0.5 - 0.5*npads)*par->GetPadPitchWidth();
1088         c.fZ = par->GetZWidth()*c.fZ; 
1089         c.fZ -= 3.*par->GetZSigma(); // PASA delay 
1090         c.fZ = sign*(z_end - c.fZ);
1091         
1092         c.fSector=sec;
1093         c.fPadRow=row;
1094         c.fTracks[0]=c.summit->fTracks[0];
1095         c.fTracks[1]=c.summit->fTracks[1];
1096         c.fTracks[2]=c.summit->fTracks[2];
1097         
1098         if (c.summit->fSignal<0) {
1099           c.fSigmaY2 *= 25.;
1100           c.fSigmaZ2 *= 4.;
1101         }
1102
1103         if (c.npeaks==0) {AddCluster(c); ncls++; ncl++;}
1104         else {
1105           new ((*fClusters)[c.idx]) AliTPCcluster(c);
1106         }
1107       }
1108     }
1109     
1110     cerr<<"sector, row, digits, clusters: "
1111         <<sec<<' '<<row<<' '<<ndigits<<' '<<ncl<<"                  \r";
1112     
1113     fDigits->Clear();
1114   
1115     delete[] bins;  
1116   }
1117 }
1118
1119 //_____________________________________________________________________________
1120 void AliTPC::ElDiff(Float_t *xyz)
1121 {
1122   //--------------------------------------------------
1123   // calculates the diffusion of a single electron
1124   //--------------------------------------------------
1125
1126   //-----------------------------------------------------------------
1127   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1128   //-----------------------------------------------------------------
1129   AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1130   Float_t driftl;
1131   
1132   Float_t z0=xyz[2];
1133
1134   driftl=z_end-TMath::Abs(xyz[2]);
1135
1136   if(driftl<0.01) driftl=0.01;
1137
1138   // check the attachment
1139
1140   driftl=TMath::Sqrt(driftl);
1141
1142   //  Float_t sig_t = driftl*diff_t;
1143   //Float_t sig_l = driftl*diff_l;
1144   Float_t sig_t = driftl*fTPCParam->GetDiffT();
1145   Float_t sig_l = driftl*fTPCParam->GetDiffL();
1146   xyz[0]=gRandom->Gaus(xyz[0],sig_t);
1147   xyz[1]=gRandom->Gaus(xyz[1],sig_t);
1148   xyz[2]=gRandom->Gaus(xyz[2],sig_l);
1149   
1150   if (TMath::Abs(xyz[2])>z_end){
1151     xyz[2]=TMath::Sign(z_end,z0);
1152   }
1153   if(xyz[2]*z0 < 0.){
1154     Float_t eps = 0.0001;
1155     xyz[2]=TMath::Sign(eps,z0);
1156   } 
1157 }
1158
1159 //_____________________________________________________________________________
1160 void AliTPC::Hits2Clusters()
1161 {
1162   //--------------------------------------------------------
1163   // TPC simple cluster generator from hits
1164   // obtained from the TPC Fast Simulator
1165   // The point errors are taken from the parametrization
1166   //--------------------------------------------------------
1167
1168   //-----------------------------------------------------------------
1169   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1170   //-----------------------------------------------------------------
1171
1172   AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1173   Float_t sigma_rphi,sigma_z,cl_rphi,cl_z;
1174   //
1175   TParticle *particle; // pointer to a given particle
1176   AliTPChit *tpcHit; // pointer to a sigle TPC hit
1177   TClonesArray *Particles; //pointer to the particle list
1178   Int_t sector,nhits;
1179   Int_t ipart;
1180   Float_t xyz[5];
1181   Float_t pl,pt,tanth,rpad,ratio;
1182   Float_t cph,sph;
1183   
1184   //---------------------------------------------------------------
1185   //  Get the access to the tracks 
1186   //---------------------------------------------------------------
1187   
1188   TTree *TH = gAlice->TreeH();
1189   Stat_t ntracks = TH->GetEntries();
1190   Particles=gAlice->Particles();
1191   
1192   //------------------------------------------------------------
1193   // Loop over all sectors (72 sectors)
1194   // Sectors 0-35 are lower sectors, 0-17 z>0, 18-35 z<0
1195   // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1196   //
1197   // First cluster for sector 0 starts at "0"
1198   //------------------------------------------------------------
1199   
1200   
1201   //
1202   int Nsectors=fDigParam->GetParam().GetNSector();
1203   for(Int_t isec=0; isec<Nsectors; isec++){
1204     //MI change
1205     fTPCParam->AdjustAngles(isec,cph,sph);
1206     
1207     //------------------------------------------------------------
1208     // Loop over tracks
1209     //------------------------------------------------------------
1210     
1211     for(Int_t track=0;track<ntracks;track++){
1212       ResetHits();
1213       TH->GetEvent(track);
1214       //
1215       //  Get number of the TPC hits and a pointer
1216       //  to the particles
1217       //
1218       nhits=fHits->GetEntriesFast();
1219       //
1220       // Loop over hits
1221       //
1222       for(Int_t hit=0;hit<nhits;hit++){
1223         tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
1224         sector=tpcHit->fSector; // sector number
1225         if(sector != isec) continue; //terminate iteration
1226         ipart=tpcHit->fTrack;
1227         particle=(TParticle*)Particles->UncheckedAt(ipart);
1228         pl=particle->Pz();
1229         pt=particle->Pt();
1230         if(pt < 1.e-9) pt=1.e-9;
1231         tanth=pl/pt;
1232         tanth = TMath::Abs(tanth);
1233         rpad=TMath::Sqrt(tpcHit->fX*tpcHit->fX + tpcHit->fY*tpcHit->fY);
1234         ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
1235         
1236         //   space-point resolutions
1237         
1238         sigma_rphi=SigmaY2(rpad,tanth,pt);
1239         sigma_z   =SigmaZ2(rpad,tanth   );
1240         
1241         //   cluster widths
1242         
1243         cl_rphi=ac_rphi-bc_rphi*rpad*tanth+cc_rphi*ratio*ratio;
1244         cl_z=ac_z-bc_z*rpad*tanth+cc_z*tanth*tanth;
1245         
1246         // temporary protection
1247         
1248         if(sigma_rphi < 0.) sigma_rphi=0.4e-3;
1249         if(sigma_z < 0.) sigma_z=0.4e-3;
1250         if(cl_rphi < 0.) cl_rphi=2.5e-3;
1251         if(cl_z < 0.) cl_z=2.5e-5;
1252         
1253         //
1254         // smearing --> rotate sectors firstly,
1255         // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1256         //
1257         Float_t xprim= tpcHit->fX*cph + tpcHit->fY*sph;
1258         Float_t yprim=-tpcHit->fX*sph + tpcHit->fY*cph;
1259         xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigma_rphi));   // y
1260         Double_t alpha=(sector < fTPCParam->GetNInnerSector()) ?
1261         fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1262         if (TMath::Abs(xyz[0]/xprim) > TMath::Tan(0.5*alpha)) xyz[0]=yprim;
1263         xyz[1]=gRandom->Gaus(tpcHit->fZ,TMath::Sqrt(sigma_z)); // z 
1264         if (TMath::Abs(xyz[1]) > 250) xyz[1]=tpcHit->fZ;
1265         xyz[2]=tpcHit->fQ+1;// q; let it be not equal to zero (Y.Belikov)
1266         xyz[3]=sigma_rphi;                                     // fSigmaY2
1267         xyz[4]=sigma_z;                                        // fSigmaZ2
1268                 
1269         // and finally add the cluster
1270         Int_t tracks[5]={tpcHit->fTrack, -1, -1, sector, tpcHit->fPadRow};
1271         AddCluster(xyz,tracks);
1272         
1273       } // end of loop over hits
1274     }   // end of loop over tracks 
1275     
1276     
1277   } // end of loop over sectors
1278   
1279   
1280 }
1281
1282
1283 void AliTPC::Hits2Digits()  
1284
1285
1286  //----------------------------------------------------
1287   // Loop over all sectors (72 sectors)
1288   // Sectors 0-35 are lower sectors, 0-17 z>0, 18-35 z<0
1289   // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1290   //----
1291   int Nsectors=fDigParam->GetParam().GetNSector();
1292   for(Int_t isec=0;isec<Nsectors;isec++)  Hits2DigitsSector(isec);
1293 }
1294
1295
1296 //_____________________________________________________________________________
1297 void AliTPC::Hits2DigitsSector(Int_t isec)
1298 {
1299   //-------------------------------------------------------------------
1300   // TPC conversion from hits to digits.
1301   //------------------------------------------------------------------- 
1302
1303   //-----------------------------------------------------------------
1304   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1305   //-----------------------------------------------------------------
1306
1307   //-------------------------------------------------------
1308   //  Get the access to the track hits
1309   //-------------------------------------------------------
1310
1311   AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1312   TTree *TH = gAlice->TreeH(); // pointer to the hits tree
1313   Stat_t ntracks = TH->GetEntries();
1314
1315   if( ntracks > 0){
1316
1317   //------------------------------------------- 
1318   //  Only if there are any tracks...
1319   //-------------------------------------------
1320
1321     
1322     // TObjArrays for three neighouring pad-rows
1323
1324     TObjArray **rowTriplet = new TObjArray* [3]; 
1325     
1326     // TObjArray-s for each pad-row
1327
1328     TObjArray **row;
1329       
1330     for (Int_t trip=0;trip<3;trip++){  
1331       rowTriplet[trip]=new TObjArray;
1332     }
1333
1334
1335     
1336       printf("*** Processing sector number %d ***\n",isec);
1337
1338       Int_t nrows =fTPCParam->GetNRow(isec);
1339
1340       row= new TObjArray* [nrows];
1341     
1342       MakeSector(isec,nrows,TH,ntracks,row);
1343
1344       //--------------------------------------------------------
1345       //   Digitize this sector, row by row
1346       //   row[i] is the pointer to the TObjArray of TVectors,
1347       //   each one containing electrons accepted on this
1348       //   row, assigned into tracks
1349       //--------------------------------------------------------
1350
1351       Int_t i;
1352
1353       for (i=0;i<nrows;i++){
1354
1355         // Triplets for i = 0 and i=1 are identical!
1356         // The same for i = nrows-1 and nrows!
1357
1358         if(i != 1 && i != nrows-1){
1359            MakeTriplet(i,rowTriplet,row);
1360          }
1361
1362             DigitizeRow(i,isec,rowTriplet);
1363
1364             fDigParam->Fill();
1365
1366             Int_t ndig=fDigParam->GetArray()->GetEntriesFast();
1367
1368             printf("*** Sector, row, digits %d %d %d ***\n",isec,i,ndig);
1369
1370             ResetDigits(); // reset digits for this row after storing them
1371              
1372        } // end of the sector digitization
1373      
1374        // delete the last triplet
1375
1376        for (i=0;i<3;i++) rowTriplet[i]->Delete();
1377           
1378        delete [] row; // delete the array of pointers to TObjArray-s
1379         
1380   } // ntracks >0
1381 } // end of Hits2Digits
1382 //_____________________________________________________________________________
1383 void AliTPC::MakeTriplet(Int_t row,
1384                          TObjArray **rowTriplet, TObjArray **prow)
1385 {
1386   //------------------------------------------------------------------
1387   //  Makes the "triplet" of the neighbouring pad-row for the
1388   //  digitization including the cross-talk between the pad-rows
1389   //------------------------------------------------------------------
1390
1391   //-----------------------------------------------------------------
1392   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1393   //-----------------------------------------------------------------
1394
1395   AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1396   Float_t gasgain = fTPCParam->GetGasGain();
1397   Int_t nTracks[3];
1398
1399   Int_t nElements,nElectrons;
1400
1401   TVector *pv;
1402   TVector *tv;
1403
1404   //-------------------------------------------------------------------
1405   // pv is an "old" track, i.e. label + triplets of (x,y,z) 
1406   // for each electron
1407   //
1408   //-------------------------------------------------------------------
1409
1410
1411   Int_t i1,i2;
1412   Int_t nel,nt;
1413
1414   if(row == 0 || row == 1){
1415
1416   // create entire triplet for the first AND the second row
1417
1418     nTracks[0] = prow[0]->GetEntries();
1419     nTracks[1] = prow[1]->GetEntries();
1420     nTracks[2] = prow[2]->GetEntries();
1421
1422     for(i1=0;i1<3;i1++){
1423       nt = nTracks[i1]; // number of tracks for this row
1424
1425       for(i2=0;i2<nt;i2++){
1426         pv = (TVector*)prow[i1]->At(i2);
1427         TVector &v1 = *pv;
1428         nElements = pv->GetNrows(); 
1429         nElectrons = (nElements-1)/3;
1430
1431         tv = new TVector(4*nElectrons+1); // create TVector for a modified track
1432         TVector &v2 = *tv;
1433         v2(0)=v1(0); //track label
1434
1435         for(nel=0;nel<nElectrons;nel++){        
1436           Int_t idx1 = nel*3;
1437           Int_t idx2 = nel*4;
1438           // Avalanche, including fluctuations
1439           Int_t aval = (Int_t) (-gasgain*TMath::Log(gRandom->Rndm()));
1440           v2(idx2+1)= v1(idx1+1);
1441           v2(idx2+2)= v1(idx1+2);
1442           v2(idx2+3)= v1(idx1+3);
1443           v2(idx2+4)= (Float_t)aval; // in number of electrons!        
1444         } // end of loop over electrons
1445         //
1446         //  Add this track to a row 
1447         //
1448
1449         rowTriplet[i1]->Add(tv); 
1450
1451
1452       } // end of loop over tracks for this row
1453
1454       prow[i1]->Delete(); // remove "old" tracks
1455       delete prow[i1]; // delete  TObjArray for this row
1456       prow[i1]=0; // set pointer to NULL
1457
1458     } // end of loop over row triplets
1459
1460
1461   }
1462   else{
1463    
1464     rowTriplet[0]->Delete(); // remove old lower row
1465
1466     nTracks[0]=rowTriplet[1]->GetEntries(); // previous middle row
1467     nTracks[1]=rowTriplet[2]->GetEntries(); // previous upper row
1468     nTracks[2]=prow[row+1]->GetEntries(); // next row
1469     
1470
1471     //------------------------------------------- 
1472     //  shift new tracks downwards
1473     //-------------------------------------------
1474
1475     for(i1=0;i1<nTracks[0];i1++){
1476       pv=(TVector*)rowTriplet[1]->At(i1);
1477       rowTriplet[0]->Add(pv); 
1478     }       
1479
1480     rowTriplet[1]->Clear(); // leave tracks on the heap!!!
1481
1482     for(i1=0;i1<nTracks[1];i1++){
1483       pv=(TVector*)rowTriplet[2]->At(i1);
1484       rowTriplet[1]->Add(pv);
1485     }
1486
1487     rowTriplet[2]->Clear(); // leave tracks on the heap!!!
1488
1489     //---------------------------------------------
1490     //  Create new upper row
1491     //---------------------------------------------
1492
1493     
1494
1495     for(i1=0;i1<nTracks[2];i1++){
1496         pv = (TVector*)prow[row+1]->At(i1);
1497         TVector &v1 = *pv;
1498         nElements = pv->GetNrows();
1499         nElectrons = (nElements-1)/3;
1500
1501         tv = new TVector(4*nElectrons+1); // create TVector for a modified track
1502         TVector &v2 = *tv;
1503         v2(0)=v1(0); //track label
1504
1505         for(nel=0;nel<nElectrons;nel++){
1506         
1507           Int_t idx1 = nel*3;
1508           Int_t idx2 = nel*4;
1509           // Avalanche, including fluctuations
1510           Int_t aval = (Int_t) 
1511             (-gasgain*TMath::Log(gRandom->Rndm()));
1512           
1513           v2(idx2+1)= v1(idx1+1); 
1514           v2(idx2+2)= v1(idx1+2);
1515           v2(idx2+3)= v1(idx1+3);
1516           v2(idx2+4)= (Float_t)aval; // in number of electrons!        
1517         } // end of loop over electrons
1518
1519           rowTriplet[2]->Add(tv);
1520      
1521     } // end of loop over tracks
1522          
1523     prow[row+1]->Delete(); // delete tracks for this row
1524     delete prow[row+1];  // delete TObjArray for this row
1525     prow[row+1]=0; // set a pointer to NULL
1526
1527   }
1528
1529 }  // end of MakeTriplet
1530 //_____________________________________________________________________________
1531 void AliTPC::ExB(Float_t *xyz)
1532 {
1533   //------------------------------------------------
1534   //  ExB at the wires and wire number calulation
1535   //------------------------------------------------
1536   
1537   //-----------------------------------------------------------------
1538   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1539   //-----------------------------------------------------------------
1540    AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1541
1542    Float_t x1=xyz[0];
1543    fTPCParam->GetWire(x1);        //calculate nearest wire position
1544    Float_t dx=xyz[0]-x1;
1545    xyz[1]+=dx*fTPCParam->GetOmegaTau();
1546
1547 } // end of ExB
1548 //_____________________________________________________________________________
1549 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rowTriplet)
1550 {
1551   //-----------------------------------------------------------
1552   // Single row digitization, coupling from the neighbouring
1553   // rows taken into account
1554   //-----------------------------------------------------------
1555
1556   //-----------------------------------------------------------------
1557   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1558   //-----------------------------------------------------------------
1559  
1560
1561   AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1562   Float_t chipgain= fTPCParam->GetChipGain();
1563   Float_t zerosup = fTPCParam->GetZeroSup();
1564   Int_t nrows =fTPCParam->GetNRow(isec);
1565   const int MAXTBKT=
1566   int((z_end+6*fTPCParam->GetZSigma())/fTPCParam->GetZWidth())+1;
1567   
1568   Int_t nTracks[3];
1569   Int_t n_of_pads[3];
1570   Int_t IndexRange[4];
1571   Int_t i1;
1572   Int_t iFlag; 
1573
1574   //
1575   // iFlag = 0 -> inner row, iFlag = 1 -> the middle one, iFlag = 2 -> the outer one
1576   //
1577
1578   nTracks[0]=rowTriplet[0]->GetEntries(); // lower row
1579   nTracks[1]=rowTriplet[1]->GetEntries(); // middle row
1580   nTracks[2]=rowTriplet[2]->GetEntries(); // upper row
1581
1582     
1583   if(irow == 0){
1584     iFlag=0;
1585     n_of_pads[0]=fTPCParam->GetNPads(isec,0);
1586     n_of_pads[1]=fTPCParam->GetNPads(isec,1);
1587   }
1588   else if(irow == nrows-1){
1589      iFlag=2;
1590      n_of_pads[1]=fTPCParam->GetNPads(isec,irow-1);
1591      n_of_pads[2]=fTPCParam->GetNPads(isec,irow);
1592   }
1593   else {
1594     iFlag=1;
1595     for(i1=0;i1<3;i1++){
1596        n_of_pads[i1]=fTPCParam->GetNPads(isec,irow-1+i1);
1597     }
1598   }
1599  
1600   //
1601   //  Integrated signal for this row
1602   //  and a single track signal
1603   // 
1604    
1605   TMatrix *m1   = new TMatrix(0,n_of_pads[iFlag]-1,0,MAXTBKT-1); // integrated
1606   TMatrix *m2   = new TMatrix(0,n_of_pads[iFlag]-1,0,MAXTBKT-1); // single
1607
1608   //
1609
1610   TMatrix &Total  = *m1;
1611
1612   //  Array of pointers to the label-signal list
1613
1614   Int_t NofDigits = n_of_pads[iFlag]*MAXTBKT; // number of digits for this row
1615
1616   Float_t  **pList = new Float_t* [NofDigits]; 
1617
1618   Int_t lp;
1619
1620   for(lp=0;lp<NofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1621
1622   //
1623   //  Straight signal and cross-talk, cross-talk is integrated over all
1624   //  tracks and added to the signal at the very end
1625   //
1626    
1627
1628   for (i1=0;i1<nTracks[iFlag];i1++){
1629
1630     m2->Zero(); // clear single track signal matrix
1631   
1632     Float_t TrackLabel = 
1633       GetSignal(rowTriplet[iFlag],i1,n_of_pads[iFlag],m2,m1,IndexRange); 
1634
1635     GetList(TrackLabel,n_of_pads[iFlag],m2,IndexRange,pList);
1636
1637   }
1638
1639   // 
1640   //  Cross talk from the neighbouring pad-rows
1641   //
1642
1643   TMatrix *m3 =  new TMatrix(0,n_of_pads[iFlag]-1,0,MAXTBKT-1); // cross-talk
1644
1645   TMatrix &Cross = *m3;
1646
1647   if(iFlag == 0){
1648
1649     // cross-talk from the outer row only (first pad row)
1650
1651     GetCrossTalk(0,rowTriplet[1],nTracks[1],n_of_pads,m3);
1652
1653   }
1654   else if(iFlag == 2){
1655
1656     // cross-talk from the inner row only (last pad row)
1657
1658     GetCrossTalk(2,rowTriplet[1],nTracks[1],n_of_pads,m3);
1659
1660   }
1661   else{
1662
1663     // cross-talk from both inner and outer rows
1664
1665     GetCrossTalk(3,rowTriplet[0],nTracks[0],n_of_pads,m3); // inner
1666     GetCrossTalk(4,rowTriplet[2],nTracks[2],n_of_pads,m3); //outer
1667   }
1668
1669   Total += Cross; // add the cross-talk
1670
1671   //
1672   //  Convert analog signal to ADC counts
1673   //
1674    
1675   Int_t tracks[3];
1676   Int_t digits[5];
1677
1678
1679   for(Int_t ip=0;ip<n_of_pads[iFlag];ip++){
1680     for(Int_t it=0;it<MAXTBKT;it++){
1681
1682       Float_t q = Total(ip,it);
1683
1684       Int_t gi =it*n_of_pads[iFlag]+ip; // global index
1685
1686       q = gRandom->Gaus(q,fTPCParam->GetNoise()); // apply noise
1687       q *= (q_el*1.e15); // convert to fC
1688       q *= chipgain; // convert to mV   
1689       q *= (adc_sat/dyn_range); // convert to ADC counts  
1690
1691       if(q <zerosup) continue; // do not fill zeros
1692       if(q > adc_sat) q = adc_sat;  // saturation
1693
1694       //
1695       //  "real" signal or electronic noise (list = -1)?
1696       //    
1697
1698       for(Int_t j1=0;j1<3;j1++){
1699         tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1700       }
1701
1702       digits[0]=isec;
1703       digits[1]=irow;
1704       digits[2]=ip;
1705       digits[3]=it;
1706       digits[4]= (Int_t)q;
1707
1708       //  Add this digit
1709
1710       AddDigit(tracks,digits);
1711     
1712     } // end of loop over time buckets
1713   }  // end of lop over pads 
1714
1715   //
1716   //  This row has been digitized, delete nonused stuff
1717   //
1718
1719   for(lp=0;lp<NofDigits;lp++){
1720     if(pList[lp]) delete [] pList[lp];
1721   }
1722   
1723   delete [] pList;
1724
1725   delete m1;
1726   delete m2;
1727   delete m3;
1728
1729 } // end of DigitizeRow
1730 //_____________________________________________________________________________
1731 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, Int_t np, TMatrix *m1, TMatrix *m2,
1732                           Int_t *IndexRange)
1733 {
1734
1735   //---------------------------------------------------------------
1736   //  Calculates 2-D signal (pad,time) for a single track,
1737   //  returns a pointer to the signal matrix and the track label 
1738   //  No digitization is performed at this level!!!
1739   //---------------------------------------------------------------
1740
1741   //-----------------------------------------------------------------
1742   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1743   //-----------------------------------------------------------------
1744
1745   TVector *tv;
1746   AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1747   AliTPCPRF2D * fPRF2D = &(fDigParam->GetPRF2D());
1748   AliTPCRF1D  * fRF    = &(fDigParam->GetRF()); 
1749   const int MAXTBKT=
1750   int((z_end+6*fTPCParam->GetZSigma())/fTPCParam->GetZWidth())+1;
1751  
1752   //to make the code faster we put parameters  to the stack
1753
1754   Float_t zwidth  = fTPCParam->GetZWidth();
1755   Float_t zwidthm1  =1./zwidth;
1756
1757   tv = (TVector*)p1->At(ntr); // pointer to a track
1758   TVector &v = *tv;
1759   
1760   Float_t label = v(0);
1761
1762   Int_t CentralPad = (np-1)/2;
1763   Int_t PadNumber;
1764   Int_t nElectrons = (tv->GetNrows()-1)/4;
1765   Float_t range=((np-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth(); // pad range
1766
1767   range -= 0.5; // dead zone, 5mm from the edge, according to H.G. Fischer
1768   
1769   Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1770
1771
1772   Float_t PadSignal[7]; // signal from a single electron
1773
1774   TMatrix &signal = *m1;
1775   TMatrix &total = *m2;
1776
1777
1778   IndexRange[0]=9999; // min pad
1779   IndexRange[1]=-1; // max pad
1780   IndexRange[2]=9999; //min time
1781   IndexRange[3]=-1; // max time
1782
1783   //
1784   //  Loop over all electrons
1785   //
1786
1787   for(Int_t nel=0; nel<nElectrons; nel++){
1788    Int_t idx=nel*4;
1789    Float_t xwire = v(idx+1);
1790    Float_t y = v(idx+2);
1791    Float_t z = v(idx+3);
1792
1793
1794    Float_t absy=TMath::Abs(y);
1795    
1796    if(absy < 0.5*fTPCParam->GetPadPitchWidth()){
1797      PadNumber=CentralPad;
1798    }
1799    else if (absy < range){
1800      PadNumber=(Int_t) ((absy-0.5*fTPCParam->GetPadPitchWidth())/fTPCParam->GetPadPitchWidth()+1.);
1801      PadNumber=(Int_t) (TMath::Sign((Float_t)PadNumber, y)+CentralPad);
1802    }
1803    else continue; // electron out of pad-range , lost at the sector edge
1804     
1805    Float_t aval = (absy<range-0.5) ? v(idx+4):v(idx+4)*IneffFactor;
1806    
1807
1808    Float_t dist = y - (Float_t)(PadNumber-CentralPad)*fTPCParam->GetPadPitchWidth();
1809    for (Int_t i=0;i<7;i++){
1810      PadSignal[i]=fPRF2D->GetPRF(-dist+(i-3)*fTPCParam->GetPadPitchWidth(),xwire)*aval;
1811      PadSignal[i] *= fTPCParam->GetPadCoupling();
1812    }
1813
1814    Int_t  LeftPad = TMath::Max(0,PadNumber-3);
1815    Int_t  RightPad = TMath::Min(np-1,PadNumber+3);
1816
1817    Int_t pmin=LeftPad-PadNumber+3; // lower index of the pad_signal vector
1818    Int_t pmax=RightPad-PadNumber+3; // upper index     
1819    
1820    Float_t z_drift = z*zwidthm1;
1821    Float_t z_offset = z_drift-(Int_t)z_drift;
1822
1823    Int_t FirstBucket = (Int_t)z_drift; // numbering starts from "0"
1824
1825
1826    // loop over time bins (4 bins is enough - 3 sigma truncated Gaussian)
1827
1828    for (Int_t i2=0;i2<4;i2++){          
1829      Int_t TrueTime = FirstBucket+i2; // current time bucket
1830      Float_t dz   = (Float_t(i2)+1.-z_offset)*zwidth; 
1831      Float_t ampl = fRF->GetRF(dz); 
1832      if( (TrueTime>MAXTBKT-1) ) break; // beyond the time range
1833      
1834      IndexRange[2]=TMath::Min(IndexRange[2],TrueTime); // min time
1835      IndexRange[3]=TMath::Max(IndexRange[3],TrueTime); // max time
1836
1837      // loop over pads, from pmin to pmax
1838      for(Int_t i3=pmin;i3<=pmax;i3++){
1839        Int_t TruePad = LeftPad+i3-pmin;
1840        IndexRange[0]=TMath::Min(IndexRange[0],TruePad); // min pad
1841        IndexRange[1]=TMath::Max(IndexRange[1],TruePad); // max pad
1842        signal(TruePad,TrueTime)+=(PadSignal[i3]*ampl); // not converted to charge!!!
1843        total(TruePad,TrueTime)+=(PadSignal[i3]*ampl); // not converted to charge!!!
1844      } // end of pads loop
1845    }  // end of time loop
1846   } // end of loop over electrons
1847
1848   return label; // returns track label when finished
1849 }
1850
1851 //_____________________________________________________________________________
1852 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *IndexRange,
1853                      Float_t **pList)
1854 {
1855   //----------------------------------------------------------------------
1856   //  Updates the list of tracks contributing to digits for a given row
1857   //----------------------------------------------------------------------
1858
1859   //-----------------------------------------------------------------
1860   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1861   //-----------------------------------------------------------------
1862
1863   TMatrix &signal = *m;
1864
1865   // lop over nonzero digits
1866
1867   for(Int_t it=IndexRange[2];it<IndexRange[3]+1;it++){
1868     for(Int_t ip=IndexRange[0];ip<IndexRange[1]+1;ip++){
1869
1870
1871         // accept only the contribution larger than 500 electrons (1/2 s_noise)
1872
1873         if(signal(ip,it)<500.) continue; 
1874
1875
1876         Int_t GlobalIndex = it*np+ip; // GlobalIndex starts from 0!
1877         
1878         if(!pList[GlobalIndex]){
1879         
1880           // 
1881           // Create new list (6 elements - 3 signals and 3 labels),
1882           //
1883
1884           pList[GlobalIndex] = new Float_t [6];
1885
1886           // set list to -1 
1887
1888           *pList[GlobalIndex] = -1.;
1889           *(pList[GlobalIndex]+1) = -1.;
1890           *(pList[GlobalIndex]+2) = -1.;
1891           *(pList[GlobalIndex]+3) = -1.;
1892           *(pList[GlobalIndex]+4) = -1.;
1893           *(pList[GlobalIndex]+5) = -1.;
1894
1895
1896           *pList[GlobalIndex] = label;
1897           *(pList[GlobalIndex]+3) = signal(ip,it);
1898         }
1899         else{
1900
1901           // check the signal magnitude
1902
1903           Float_t highest = *(pList[GlobalIndex]+3);
1904           Float_t middle = *(pList[GlobalIndex]+4);
1905           Float_t lowest = *(pList[GlobalIndex]+5);
1906
1907           //
1908           //  compare the new signal with already existing list
1909           //
1910
1911           if(signal(ip,it)<lowest) continue; // neglect this track
1912
1913           //
1914
1915           if (signal(ip,it)>highest){
1916             *(pList[GlobalIndex]+5) = middle;
1917             *(pList[GlobalIndex]+4) = highest;
1918             *(pList[GlobalIndex]+3) = signal(ip,it);
1919
1920             *(pList[GlobalIndex]+2) = *(pList[GlobalIndex]+1);
1921             *(pList[GlobalIndex]+1) = *pList[GlobalIndex];
1922             *pList[GlobalIndex] = label;
1923           }
1924           else if (signal(ip,it)>middle){
1925             *(pList[GlobalIndex]+5) = middle;
1926             *(pList[GlobalIndex]+4) = signal(ip,it);
1927
1928             *(pList[GlobalIndex]+2) = *(pList[GlobalIndex]+1);
1929             *(pList[GlobalIndex]+1) = label;
1930           }
1931           else{
1932             *(pList[GlobalIndex]+5) = signal(ip,it);
1933             *(pList[GlobalIndex]+2) = label;
1934           }
1935         }
1936
1937     } // end of loop over pads
1938   } // end of loop over time bins
1939
1940
1941
1942
1943 }//end of GetList
1944 //___________________________________________________________________
1945 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1946                         Stat_t ntracks,TObjArray **row)
1947 {
1948
1949   //-----------------------------------------------------------------
1950   // Prepares the sector digitization, creates the vectors of
1951   // tracks for each row of this sector. The track vector
1952   // contains the track label and the position of electrons.
1953   //-----------------------------------------------------------------
1954
1955   //-----------------------------------------------------------------
1956   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1957   //-----------------------------------------------------------------
1958
1959   AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1960   Int_t i;
1961   Float_t xyz[3]; 
1962
1963   AliTPChit *tpcHit; // pointer to a sigle TPC hit    
1964  
1965   //----------------------------------------------
1966   // Create TObjArray-s, one for each row,
1967   // each TObjArray will store the TVectors
1968   // of electrons, one TVector per each track.
1969   //---------------------------------------------- 
1970     
1971   for(i=0; i<nrows; i++){
1972     row[i] = new TObjArray;
1973   }
1974   Int_t *n_of_electrons = new Int_t [nrows]; // electron counter for each row
1975   TVector **tr = new TVector* [nrows]; //pointers to the track vectors
1976
1977   //--------------------------------------------------------------------
1978   //  Loop over tracks, the "track" contains the full history
1979   //--------------------------------------------------------------------
1980
1981   Int_t previousTrack,currentTrack;
1982   previousTrack = -1; // nothing to store so far!
1983
1984   for(Int_t track=0;track<ntracks;track++){
1985
1986     ResetHits();
1987
1988     TH->GetEvent(track); // get next track
1989     Int_t nhits = fHits->GetEntriesFast(); // get number of hits for this track
1990
1991     if(nhits == 0) continue; // no hits in the TPC for this track
1992
1993     //--------------------------------------------------------------
1994     //  Loop over hits
1995     //--------------------------------------------------------------
1996
1997     for(Int_t hit=0;hit<nhits;hit++){
1998
1999       tpcHit = (AliTPChit*)fHits->UncheckedAt(hit); // get a pointer to a hit
2000       
2001       Int_t sector=tpcHit->fSector; // sector number
2002       if(sector != isec) continue; //terminate iteration
2003
2004         currentTrack = tpcHit->fTrack; // track number
2005         if(currentTrack != previousTrack){
2006                           
2007            // store already filled fTrack
2008               
2009            for(i=0;i<nrows;i++){
2010              if(previousTrack != -1){
2011                if(n_of_electrons[i]>0){
2012                  TVector &v = *tr[i];
2013                  v(0) = previousTrack;
2014                  tr[i]->ResizeTo(3*n_of_electrons[i]+1); // shrink if necessary
2015                  row[i]->Add(tr[i]);                     
2016                }
2017                else{
2018                  delete tr[i]; // delete empty TVector
2019                  tr[i]=0;
2020                }
2021              }
2022
2023              n_of_electrons[i]=0;
2024              tr[i] = new TVector(361); // TVectors for the next fTrack
2025
2026            } // end of loop over rows
2027                
2028            previousTrack=currentTrack; // update track label 
2029         }
2030            
2031         Int_t QI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2032
2033        //---------------------------------------------------
2034        //  Calculate the electron attachment probability
2035        //---------------------------------------------------
2036
2037         Float_t time = 1.e6*(z_end-TMath::Abs(tpcHit->fZ))/fTPCParam->GetDriftV(); 
2038         // in microseconds!     
2039         Float_t AttProb = fTPCParam->GetAttCoef()*
2040           fTPCParam->GetOxyCont()*time; //  fraction! 
2041    
2042         //-----------------------------------------------
2043         //  Loop over electrons
2044         //-----------------------------------------------
2045
2046         for(Int_t nel=0;nel<QI;nel++){
2047           // skip if electron lost due to the attachment
2048           if((gRandom->Rndm(0)) < AttProb) continue; // electron lost!
2049           xyz[0]=tpcHit->fX;
2050           xyz[1]=tpcHit->fY;
2051           xyz[2]=tpcHit->fZ;
2052           ElDiff(xyz); // Appply the diffusion
2053           Int_t row_number;
2054           fTPCParam->XYZtoCRXYZ(xyz,isec,row_number,3);
2055
2056           //transform position to local coordinates
2057           //option 3 means that we calculate x position relative to 
2058           //nearest pad row 
2059
2060           if ((row_number<0)||row_number>=fTPCParam->GetNRow(isec)) continue;
2061           ExB(xyz); // ExB effect at the sense wires
2062           n_of_electrons[row_number]++;   
2063           //----------------------------------
2064           // Expand vector if necessary
2065           //----------------------------------
2066           if(n_of_electrons[row_number]>120){
2067             Int_t range = tr[row_number]->GetNrows();
2068             if(n_of_electrons[row_number] > (range-1)/3){
2069               tr[row_number]->ResizeTo(range+150); // Add 50 electrons
2070             }
2071           }
2072           
2073           TVector &v = *tr[row_number];
2074           Int_t idx = 3*n_of_electrons[row_number]-2;
2075
2076           v(idx)=  xyz[0];   // X
2077           v(idx+1)=xyz[1];   // Y (along the pad-row)
2078           v(idx+2)=xyz[2];   // Z
2079             
2080         } // end of loop over electrons
2081         
2082       } // end of loop over hits
2083     } // end of loop over tracks
2084
2085     //
2086     //   store remaining track (the last one) if not empty
2087     //
2088
2089      for(i=0;i<nrows;i++){
2090        if(n_of_electrons[i]>0){
2091           TVector &v = *tr[i];
2092           v(0) = previousTrack;
2093           tr[i]->ResizeTo(3*n_of_electrons[i]+1); // shrink if necessary
2094           row[i]->Add(tr[i]);  
2095         }
2096         else{
2097           delete tr[i];
2098           tr[i]=0;
2099         }  
2100       }  
2101
2102           delete [] tr;
2103           delete [] n_of_electrons; 
2104
2105 } // end of MakeSector
2106 //_____________________________________________________________________________
2107 void AliTPC::GetCrossTalk (Int_t iFlag,TObjArray *p,Int_t ntracks,Int_t *npads,
2108                            TMatrix *m)
2109 {
2110
2111   //-------------------------------------------------------------
2112   // Calculates the cross-talk from one row (inner or outer)
2113   //-------------------------------------------------------------
2114
2115   //-----------------------------------------------------------------
2116   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2117   //-----------------------------------------------------------------
2118
2119   //
2120   // iFlag=2 & 3 --> cross-talk from the inner row
2121   // iFlag=0 & 4 --> cross-talk from the outer row
2122   //
2123   AliTPCParam * fTPCParam = &(fDigParam->GetParam());
2124   AliTPCPRF2D * fPRF2D = &(fDigParam->GetPRF2D());
2125   AliTPCRF1D  * fRF    = &(fDigParam->GetRF()); 
2126   const int MAXTBKT=
2127   int((z_end+6*fTPCParam->GetZSigma())/fTPCParam->GetZWidth())+1;
2128  
2129   //to make code faster
2130
2131   Float_t zwidth  = fTPCParam->GetZWidth();
2132   Float_t zwidthm1  =1/fTPCParam->GetZWidth();
2133
2134  Int_t PadNumber;
2135  Float_t xwire;
2136
2137  Int_t nPadsSignal; // for this pads the signal is calculated
2138  Float_t range;     // sense wire range
2139  Int_t nPadsDiff;
2140
2141  Float_t IneffFactor=0.5; // gain inefficiency close to the edges
2142
2143  if(iFlag == 0){   
2144    // 1-->0
2145    nPadsSignal = npads[1];
2146    range = ((npads[1]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();  
2147    nPadsDiff = (npads[1]-npads[0])/2;
2148  }  
2149  else if (iFlag == 2){
2150    // 1-->2
2151    nPadsSignal = npads[2];
2152    range = ((npads[1]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2153    nPadsDiff = 0;
2154  }
2155  else if (iFlag == 3){
2156    // 0-->1
2157    nPadsSignal = npads[1];
2158    range = ((npads[0]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2159    nPadsDiff = 0;
2160  }
2161  else{
2162    // 2-->1
2163    nPadsSignal = npads[2];
2164    range = ((npads[2]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2165    nPadsDiff = (npads[2]-npads[1])/2;
2166  }
2167
2168  range-=0.5; // dead zone close to the edges
2169
2170  TVector *tv; 
2171  TMatrix &signal = *m;
2172
2173  Int_t CentralPad = (nPadsSignal-1)/2;
2174  Float_t PadSignal[7]; // signal from a single electron
2175  // Loop over tracks
2176  for(Int_t nt=0;nt<ntracks;nt++){
2177    tv=(TVector*)p->At(nt); // pointer to a track
2178    TVector &v = *tv;
2179    Int_t nElectrons = (tv->GetNrows()-1)/4;
2180    // Loop over electrons
2181    for(Int_t nel=0; nel<nElectrons; nel++){
2182      Int_t idx=nel*4;
2183      xwire=v(idx+1);
2184  
2185      if (iFlag==0) xwire+=fTPCParam->GetPadPitchLength();
2186      if (iFlag==2)  xwire-=fTPCParam->GetPadPitchLength();
2187      if (iFlag==3)  xwire-=fTPCParam->GetPadPitchLength();
2188      if (iFlag==4)  xwire+=fTPCParam->GetPadPitchLength();  
2189    
2190      //  electron acceptance for the cross-talk (only the closest wire)  
2191
2192      Float_t dxMax = fTPCParam->GetPadPitchLength()*0.5+fTPCParam->GetWWPitch();
2193      if(TMath::Abs(xwire)>dxMax) continue;
2194
2195      Float_t y = v(idx+2);
2196      Float_t z = v(idx+3);
2197      Float_t absy=TMath::Abs(y);
2198
2199      if(absy < 0.5*fTPCParam->GetPadPitchWidth()){
2200        PadNumber=CentralPad;
2201      }
2202      else if (absy < range){
2203        PadNumber=(Int_t) ((absy-0.5*fTPCParam->GetPadPitchWidth())/fTPCParam->GetPadPitchWidth() +1.);
2204        PadNumber=(Int_t) (TMath::Sign((Float_t)PadNumber, y)+CentralPad);
2205      }
2206      else continue; // electron out of sense wire range, lost at the sector edge
2207
2208      Float_t aval = (absy<range-0.5) ? v(idx+4):v(idx+4)*IneffFactor;
2209
2210      Float_t dist = y - (Float_t)(PadNumber-CentralPad)*fTPCParam->GetPadPitchWidth();
2211        
2212      for (Int_t i=0;i<7;i++){
2213        PadSignal[i]=fPRF2D->GetPRF(-dist+(i-3)*fTPCParam->GetPadPitchWidth(),xwire)*aval;
2214
2215        PadSignal[i] *= fTPCParam->GetPadCoupling();
2216      }
2217      // real pad range
2218
2219      Int_t  LeftPad = TMath::Max(0,PadNumber-3);
2220      Int_t  RightPad = TMath::Min(nPadsSignal-1,PadNumber+3);
2221
2222      Int_t pmin=LeftPad-PadNumber+3; // lower index of the pad_signal vector
2223      Int_t pmax=RightPad-PadNumber+3; // upper index  
2224
2225
2226      Float_t z_drift = z*zwidthm1;
2227      Float_t z_offset = z_drift-(Int_t)z_drift;
2228
2229      Int_t FirstBucket = (Int_t)z_drift; // numbering starts from "0"
2230
2231      for (Int_t i2=0;i2<4;i2++){     
2232        Int_t TrueTime = FirstBucket+i2; // current time bucket
2233        Float_t dz   = (Float_t(i2)+1.- z_offset)*zwidth; 
2234        Float_t ampl = fRF->GetRF(dz); 
2235        if((TrueTime>MAXTBKT-1)) break; // beyond the time range
2236
2237
2238        // loop over pads, from pmin to pmax
2239
2240        for(Int_t i3=pmin;i3<pmax+1;i3++){
2241          Int_t TruePad = LeftPad+i3-pmin;
2242
2243          if(TruePad<nPadsDiff || TruePad > nPadsSignal-nPadsDiff-1) continue;
2244
2245          TruePad -= nPadsDiff;
2246          signal(TruePad,TrueTime)+=(PadSignal[i3]*ampl); // not converted to charge!
2247
2248        } // end of loop over pads
2249      } // end of loop over time bins
2250
2251    } // end of loop over electrons 
2252
2253  } // end of loop over tracks
2254
2255 } // end of GetCrossTalk
2256
2257
2258
2259 //_____________________________________________________________________________
2260 void AliTPC::Init()
2261 {
2262   //
2263   // Initialise TPC detector after definition of geometry
2264   //
2265   Int_t i;
2266   //
2267   printf("\n");
2268   for(i=0;i<35;i++) printf("*");
2269   printf(" TPC_INIT ");
2270   for(i=0;i<35;i++) printf("*");
2271   printf("\n");
2272   //
2273   for(i=0;i<80;i++) printf("*");
2274   printf("\n");
2275 }
2276
2277 //_____________________________________________________________________________
2278 void AliTPC::MakeBranch(Option_t* option)
2279 {
2280   //
2281   // Create Tree branches for the TPC.
2282   //
2283   Int_t buffersize = 4000;
2284   char branchname[10];
2285   sprintf(branchname,"%s",GetName());
2286
2287   AliDetector::MakeBranch(option);
2288
2289   char *D = strstr(option,"D");
2290
2291   if (fDigits   && gAlice->TreeD() && D) {
2292     gAlice->TreeD()->Branch(branchname,&fDigits, buffersize);
2293     printf("Making Branch %s for digits\n",branchname);
2294   }     
2295
2296   char *R = strstr(option,"R");
2297
2298   if (fClusters && gAlice->TreeR() && R) {
2299     gAlice->TreeR()->Branch(branchname,&fClusters, buffersize);
2300     printf("Making Branch %s for Clusters\n",branchname);
2301   }     
2302 }
2303  
2304 //_____________________________________________________________________________
2305 void AliTPC::ResetDigits()
2306 {
2307   //
2308   // Reset number of digits and the digits array for this detector
2309   // reset clusters
2310   //
2311   fNdigits   = 0;
2312   //  if (fDigits)   fDigits->Clear();
2313   if (fDigParam->GetArray()!=0)  fDigParam->GetArray()->Clear();
2314   fNclusters = 0;
2315   if (fClusters) fClusters->Clear();
2316 }
2317
2318 //_____________________________________________________________________________
2319 void AliTPC::SetSecAL(Int_t sec)
2320 {
2321   //---------------------------------------------------
2322   // Activate/deactivate selection for lower sectors
2323   //---------------------------------------------------
2324
2325   //-----------------------------------------------------------------
2326   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2327   //-----------------------------------------------------------------
2328
2329   fSecAL = sec;
2330 }
2331
2332 //_____________________________________________________________________________
2333 void AliTPC::SetSecAU(Int_t sec)
2334 {
2335   //----------------------------------------------------
2336   // Activate/deactivate selection for upper sectors
2337   //---------------------------------------------------
2338
2339   //-----------------------------------------------------------------
2340   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2341   //-----------------------------------------------------------------
2342
2343   fSecAU = sec;
2344 }
2345
2346 //_____________________________________________________________________________
2347 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2348 {
2349   //----------------------------------------
2350   // Select active lower sectors
2351   //----------------------------------------
2352
2353   //-----------------------------------------------------------------
2354   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2355   //-----------------------------------------------------------------
2356
2357   fSecLows[0] = s1;
2358   fSecLows[1] = s2;
2359   fSecLows[2] = s3;
2360   fSecLows[3] = s4;
2361   fSecLows[4] = s5;
2362   fSecLows[5] = s6;
2363 }
2364
2365 //_____________________________________________________________________________
2366 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2367                        Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10, 
2368                        Int_t s11 , Int_t s12)
2369 {
2370   //--------------------------------
2371   // Select active upper sectors
2372   //--------------------------------
2373
2374   //-----------------------------------------------------------------
2375   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2376   //-----------------------------------------------------------------
2377
2378   fSecUps[0] = s1;
2379   fSecUps[1] = s2;
2380   fSecUps[2] = s3;
2381   fSecUps[3] = s4;
2382   fSecUps[4] = s5;
2383   fSecUps[5] = s6;
2384   fSecUps[6] = s7;
2385   fSecUps[7] = s8;
2386   fSecUps[8] = s9;
2387   fSecUps[9] = s10;
2388   fSecUps[10] = s11;
2389   fSecUps[11] = s12;
2390 }
2391
2392 //_____________________________________________________________________________
2393 void AliTPC::SetSens(Int_t sens)
2394 {
2395
2396   //-------------------------------------------------------------
2397   // Activates/deactivates the sensitive strips at the center of
2398   // the pad row -- this is for the space-point resolution calculations
2399   //-------------------------------------------------------------
2400
2401   //-----------------------------------------------------------------
2402   // Origin: Marek Kowalski  IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2403   //-----------------------------------------------------------------
2404
2405   fSens = sens;
2406 }
2407  
2408 void AliTPC::SetSide(Float_t side)
2409 {
2410   fSide = side;
2411  
2412 }
2413 //____________________________________________________________________________
2414 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2415                            Float_t p2,Float_t p3)
2416 {
2417
2418  fNoComp = nc;
2419  
2420  fMixtComp[0]=c1;
2421  fMixtComp[1]=c2;
2422  fMixtComp[2]=c3;
2423
2424  fMixtProp[0]=p1;
2425  fMixtProp[1]=p2;
2426  fMixtProp[2]=p3; 
2427  
2428  
2429 }
2430 //_____________________________________________________________________________
2431 void AliTPC::Streamer(TBuffer &R__b)
2432 {
2433   //
2434   // Stream an object of class AliTPC.
2435   //
2436    if (R__b.IsReading()) {
2437       Version_t R__v = R__b.ReadVersion(); if (R__v) { }
2438       AliDetector::Streamer(R__b);
2439       if (R__v < 2) return;
2440       R__b >> fNsectors;
2441       R__b >> fNclusters;
2442       R__b >> fNtracks;
2443       fDigitsIndex   = new Int_t[fNsectors+1];
2444    } else {
2445       R__b.WriteVersion(AliTPC::IsA());
2446       AliDetector::Streamer(R__b);
2447       R__b << fNsectors;
2448       R__b << fNclusters;
2449       R__b << fNtracks;
2450    }
2451 }
2452   
2453 ClassImp(AliTPCcluster)
2454  
2455 //_____________________________________________________________________________
2456 AliTPCcluster::AliTPCcluster(Float_t *hits, Int_t *lab)
2457 {
2458   //
2459   // Creates a simulated cluster for the TPC
2460   //
2461   fTracks[0]  = lab[0];
2462   fTracks[1]  = lab[1];
2463   fTracks[2]  = lab[2];
2464   fSector     = lab[3];
2465   fPadRow     = lab[4];
2466   fY          = hits[0];
2467   fZ          = hits[1];
2468   fQ          = hits[2];
2469   fSigmaY2    = hits[3];
2470   fSigmaZ2    = hits[4];
2471 }
2472  
2473 //_____________________________________________________________________________
2474 void AliTPCcluster::GetXYZ(Float_t *x, const AliTPCParam *par) const 
2475 {
2476   //
2477   // Transformation from local to global coordinate system
2478   //
2479   x[0]=par->GetPadRowRadii(fSector,fPadRow);
2480   x[1]=fY;
2481   x[2]=fZ;
2482   par->CRXYZtoXYZ(x,fSector,fPadRow,1);
2483   x[2]=fZ;
2484 }
2485  
2486 //_____________________________________________________________________________
2487 Int_t AliTPCcluster::Compare(TObject * o)
2488 {
2489   //
2490   // compare two clusters according y coordinata
2491   //
2492   AliTPCcluster *cl= (AliTPCcluster *)o;
2493   if (fY<cl->fY) return -1;
2494   if (fY==cl->fY) return 0;
2495   return 1;  
2496 }
2497
2498 Bool_t AliTPCcluster::IsSortable() const
2499 {
2500   //
2501   //make AliTPCcluster sortabale
2502   //
2503   return kTRUE; 
2504 }
2505
2506
2507
2508 ClassImp(AliTPCdigit)
2509  
2510 //_____________________________________________________________________________
2511 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2512   AliDigit(tracks)
2513 {
2514   //
2515   // Creates a TPC digit object
2516   //
2517   fSector     = digits[0];
2518   fPadRow     = digits[1];
2519   fPad        = digits[2];
2520   fTime       = digits[3];
2521   fSignal     = digits[4];
2522 }
2523
2524  
2525 ClassImp(AliTPChit)
2526  
2527 //_____________________________________________________________________________
2528 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2529 AliHit(shunt,track)
2530 {
2531   //
2532   // Creates a TPC hit object
2533   //
2534   fSector     = vol[0];
2535   fPadRow     = vol[1];
2536   fX          = hits[0];
2537   fY          = hits[1];
2538   fZ          = hits[2];
2539   fQ          = hits[3];
2540 }
2541  
2542  
2543 ClassImp(AliTPCtrack)
2544  
2545 //_____________________________________________________________________________
2546 AliTPCtrack::AliTPCtrack(Float_t *hits)
2547 {
2548   //
2549   // Default creator for a TPC reconstructed track object
2550   //
2551   fX=hits[0]; // This is dummy code !
2552 }
2553
2554 AliTPCtrack::AliTPCtrack(const AliTPCcluster *c,const TVector& xx,
2555                          const TMatrix& CC, Double_t xref, Double_t alpha):
2556   x(xx),C(CC),fClusters(200)
2557 {
2558   //-----------------------------------------------------------------
2559   // This is the main track constructor.
2560   //
2561   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2562   //-----------------------------------------------------------------
2563   fX=xref;
2564   fAlpha=alpha;
2565   fChi2=0.;
2566   fClusters.AddLast((AliTPCcluster*)(c));
2567 }
2568
2569 //_____________________________________________________________________________
2570 AliTPCtrack::AliTPCtrack(const AliTPCtrack& t) : x(t.x), C(t.C),
2571   fClusters(t.fClusters.GetEntriesFast()) 
2572 {
2573   //-----------------------------------------------------------------
2574   // This is a track copy constructor.
2575   //
2576   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2577   //-----------------------------------------------------------------
2578   fX=t.fX;
2579   fChi2=t.fChi2;
2580   fAlpha=t.fAlpha;
2581   int n=t.fClusters.GetEntriesFast();
2582   for (int i=0; i<n; i++) fClusters.AddLast(t.fClusters.UncheckedAt(i));
2583 }
2584
2585 //_____________________________________________________________________________
2586 Int_t AliTPCtrack::Compare(TObject *o) {
2587   //-----------------------------------------------------------------
2588   // This function compares tracks according to their curvature.
2589   //
2590   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2591   //-----------------------------------------------------------------
2592   AliTPCtrack *t=(AliTPCtrack*)o;
2593   Double_t co=t->GetSigmaY2();
2594   Double_t c =GetSigmaY2();
2595   if (c>co) return 1;
2596   else if (c<co) return -1;
2597   return 0;
2598 }
2599
2600 //_____________________________________________________________________________
2601 int AliTPCtrack::PropagateTo(Double_t xk,Double_t x0,Double_t rho,Double_t pm)
2602 {
2603   //-----------------------------------------------------------------
2604   // This function propagates a track to a reference plane x=xk.
2605   //
2606   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2607   //-----------------------------------------------------------------
2608   if (TMath::Abs(x(2)*xk - x(3)) >= 0.999) {
2609     if (*this>4) cerr<<*this<<" AliTPCtrack warning: Propagation failed !\n";
2610     return 0;
2611   }
2612
2613   Double_t x1=fX, x2=x1+0.5*(xk-x1), dx=x2-x1, y1=x(0), z1=x(1);
2614   Double_t c1=x(2)*x1 - x(3), r1=sqrt(1.- c1*c1);
2615   Double_t c2=x(2)*x2 - x(3), r2=sqrt(1.- c2*c2);
2616   
2617   x(0) += dx*(c1+c2)/(r1+r2);
2618   x(1) += dx*(c1+c2)/(c1*r2 + c2*r1)*x(4);
2619
2620   TMatrix F(5,5); F.UnitMatrix();
2621   Double_t rr=r1+r2, cc=c1+c2, xx=x1+x2;
2622   F(0,2)= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
2623   F(0,3)=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
2624   Double_t cr=c1*r2+c2*r1;
2625   F(1,2)= dx*x(4)*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
2626   F(1,3)=-dx*x(4)*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
2627   F(1,4)= dx*cc/cr; 
2628   TMatrix tmp(F,TMatrix::kMult,C);
2629   C.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2630   
2631   fX=x2;
2632   
2633   //Multiple scattering******************
2634   Double_t ey=x(2)*fX - x(3);
2635   Double_t ex=sqrt(1-ey*ey);
2636   Double_t ez=x(4);
2637   TMatrix Q(5,5); Q=0.;
2638   Q(2,2)=ez*ez+ey*ey;   Q(2,3)=-ex*ey;       Q(2,4)=-ex*ez;
2639   Q(3,2)=Q(2,3);        Q(3,3)= ez*ez+ex*ex; Q(3,4)=-ey*ez;
2640   Q(4,2)=Q(2,4);        Q(4,3)= Q(3,4);      Q(4,4)=1.;
2641   
2642   F=0;
2643   F(2,2)=-x(2)*ex;          F(2,3)=-x(2)*ey;
2644   F(3,2)=-ex*(x(2)*fX-ey);  F(3,3)=-(1.+ x(2)*fX*ey - ey*ey);
2645   F(4,2)=-ez*ex;            F(4,3)=-ez*ey;           F(4,4)=1.;
2646   
2647   tmp.Mult(F,Q);
2648   Q.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2649   
2650   Double_t p2=GetPt()*GetPt()*(1.+x(4)*x(4));
2651   Double_t beta2=p2/(p2 + pm*pm);
2652   Double_t d=sqrt((x1-fX)*(x1-fX)+(y1-x(0))*(y1-x(0))+(z1-x(1))*(z1-x(1)));
2653   d*=2.;
2654   Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho;
2655   Q*=theta2;
2656   C+=Q;
2657   
2658   //Energy losses************************
2659   Double_t dE=0.153e-3/beta2*(log(5940*beta2/(1-beta2)) - beta2)*d*rho;
2660   if (x1 < x2) dE=-dE;
2661   x(2)*=(1.- sqrt(p2+pm*pm)/p2*dE);
2662   //x(3)*=(1.- sqrt(p2+pm*pm)/p2*dE);
2663   
2664   x1=fX; x2=xk; y1=x(0); z1=x(1);
2665   c1=x(2)*x1 - x(3); r1=sqrt(1.- c1*c1);
2666   c2=x(2)*x2 - x(3); r2=sqrt(1.- c2*c2);
2667   
2668   x(0) += dx*(c1+c2)/(r1+r2);
2669   x(1) += dx*(c1+c2)/(c1*r2 + c2*r1)*x(4);
2670   
2671   F.UnitMatrix();
2672   rr=r1+r2; cc=c1+c2; xx=x1+x2;
2673   F(0,2)= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
2674   F(0,3)=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
2675   cr=c1*r2+c2*r1;
2676   F(1,2)= dx*x(4)*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
2677   F(1,3)=-dx*x(4)*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
2678   F(1,4)= dx*cc/cr; 
2679   tmp.Mult(F,C);
2680   C.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2681   
2682   fX=x2;
2683   
2684   return 1;
2685 }
2686
2687 //_____________________________________________________________________________
2688 void AliTPCtrack::PropagateToVertex(Double_t x0,Double_t rho,Double_t pm) 
2689 {
2690   //-----------------------------------------------------------------
2691   // This function propagates tracks to the "vertex".
2692   //
2693   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2694   //-----------------------------------------------------------------
2695   Double_t c=x(2)*fX - x(3);
2696   Double_t tgf=-x(3)/(x(2)*x(0) + sqrt(1-c*c));
2697   Double_t snf=tgf/sqrt(1.+ tgf*tgf);
2698   Double_t xv=(x(3)+snf)/x(2);
2699   PropagateTo(xv,x0,rho,pm);
2700 }
2701
2702 //_____________________________________________________________________________
2703 void AliTPCtrack::Update(const AliTPCcluster *c, Double_t chisq)
2704 {
2705   //-----------------------------------------------------------------
2706   // This function associates a clusters with this track.
2707   //
2708   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2709   //-----------------------------------------------------------------
2710   TMatrix H(2,5); H.UnitMatrix();
2711   TMatrix Ht(TMatrix::kTransposed,H);
2712   TVector m(2);   m(0)=c->fY; m(1)=c->fZ; 
2713   TMatrix V(2,2); V(0,0)=c->fSigmaY2; V(0,1)=0.; V(1,0)=0.; V(1,1)=c->fSigmaZ2;
2714
2715   TMatrix tmp(H,TMatrix::kMult,C);
2716   TMatrix R(tmp,TMatrix::kMult,Ht); R+=V;
2717   
2718   Double_t det=(Double_t)R(0,0)*R(1,1) - (Double_t)R(0,1)*R(1,0);
2719   R(0,1)=R(0,0); R(0,0)=R(1,1); R(1,1)=R(0,1); 
2720   R(1,0)*=-1; R(0,1)=R(1,0);
2721   R*=1./det;
2722   
2723   //R.Invert();
2724   
2725   TMatrix K(C,TMatrix::kMult,Ht); K*=R;
2726   
2727   TVector savex=x;
2728   x*=H; x-=m; x*=-1; x*=K; x+=savex;
2729   if (TMath::Abs(x(2)*fX-x(3)) >= 0.999) {
2730     if (*this>4) cerr<<*this<<" AliTPCtrack warning: Filtering failed !\n";
2731     x=savex;
2732     return;
2733   }
2734   
2735   TMatrix saveC=C;
2736   C.Mult(K,tmp); C-=saveC; C*=-1;
2737   
2738   fClusters.AddLast((AliTPCcluster*)c);
2739   fChi2 += chisq;
2740 }
2741
2742 //_____________________________________________________________________________
2743 int AliTPCtrack::Rotate(Double_t alpha)
2744 {
2745   //-----------------------------------------------------------------
2746   // This function rotates this track.
2747   //
2748   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2749   //-----------------------------------------------------------------
2750   fAlpha += alpha;
2751   
2752   Double_t x1=fX, y1=x(0);
2753   Double_t ca=cos(alpha), sa=sin(alpha);
2754   Double_t r1=x(2)*fX - x(3);
2755   
2756   fX = x1*ca + y1*sa;
2757   x(0)=-x1*sa + y1*ca;
2758   x(3)=x(3)*ca + (x(2)*y1 + sqrt(1.- r1*r1))*sa;
2759   
2760   Double_t r2=x(2)*fX - x(3);
2761   if (TMath::Abs(r2) >= 0.999) {
2762     if (*this>4) cerr<<*this<<" AliTPCtrack warning: Rotation failed !\n";
2763     return 0;
2764   }
2765   
2766   Double_t y0=x(0) + sqrt(1.- r2*r2)/x(2);
2767   if ((x(0)-y0)*x(2) >= 0.) {
2768     if (*this>4) cerr<<*this<<" AliTPCtrack warning: Rotation failed !!!\n";
2769     return 0;
2770   }
2771   
2772   TMatrix F(5,5); F.UnitMatrix();
2773   F(0,0)=ca;
2774   F(3,0)=x(2)*sa; 
2775   F(3,2)=(y1 - r1*x1/sqrt(1.- r1*r1))*sa; 
2776   F(3,3)= ca + sa*r1/sqrt(1.- r1*r1);
2777   TMatrix tmp(F,TMatrix::kMult,C); 
2778   // Double_t dy2=C(0,0);
2779   C.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2780   // C(0,0)+=dy2*sa*sa*r1*r1/(1.- r1*r1);
2781   // C(1,1)+=dy2*sa*sa*x(4)*x(4)/(1.- r1*r1);
2782   
2783   return 1;
2784 }
2785
2786 //_____________________________________________________________________________
2787 void AliTPCtrack::UseClusters() const 
2788 {
2789   //-----------------------------------------------------------------
2790   // This function marks clusters associated with this track.
2791   //
2792   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2793   //-----------------------------------------------------------------
2794   int num_of_clusters=fClusters.GetEntriesFast();
2795   for (int i=0; i<num_of_clusters; i++) {
2796     //if (i<=14) continue;
2797     AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(i);
2798     c->Use();   
2799   }
2800 }
2801
2802 //_____________________________________________________________________________
2803 Double_t AliTPCtrack::GetPredictedChi2(const AliTPCcluster *c) const 
2804 {
2805   //-----------------------------------------------------------------
2806   // This function calculates a predicted chi2 increment.
2807   //
2808   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2809   //-----------------------------------------------------------------
2810   TMatrix H(2,5); H.UnitMatrix();
2811   TVector m(2);   m(0)=c->fY; m(1)=c->fZ; 
2812   TMatrix V(2,2); V(0,0)=c->fSigmaY2; V(0,1)=0.; V(1,0)=0.; V(1,1)=c->fSigmaZ2;
2813   TVector res=x;  res*=H; res-=m; //res*=-1; 
2814   TMatrix tmp(H,TMatrix::kMult,C);
2815   TMatrix R(tmp,TMatrix::kMult,TMatrix(TMatrix::kTransposed,H)); R+=V;
2816   
2817   Double_t det=(Double_t)R(0,0)*R(1,1) - (Double_t)R(0,1)*R(1,0);
2818   if (TMath::Abs(det) < 1.e-10) {
2819     if (*this>4) cerr<<*this<<" AliTPCtrack warning: Singular matrix !\n";
2820     return 1e10;
2821   }
2822   R(0,1)=R(0,0); R(0,0)=R(1,1); R(1,1)=R(0,1); 
2823   R(1,0)*=-1; R(0,1)=R(1,0);
2824   R*=1./det;
2825   
2826   //R.Invert();
2827   
2828   TVector r=res;
2829   res*=R;
2830   return r*res;
2831 }
2832
2833 //_____________________________________________________________________________
2834 struct S { int lab; int max; };
2835 int AliTPCtrack::GetLabel(int nrows) const 
2836 {
2837   //-----------------------------------------------------------------
2838   // This function returns the track label. If label<0, this track is fake.
2839   //
2840   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2841   //-----------------------------------------------------------------
2842   int num_of_clusters=fClusters.GetEntriesFast();
2843   S *s=new S[num_of_clusters];
2844   int i;
2845   for (i=0; i<num_of_clusters; i++) s[i].lab=s[i].max=0;
2846   
2847   int lab=123456789;
2848   for (i=0; i<num_of_clusters; i++) {
2849     AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(i);
2850     lab=TMath::Abs(c->fTracks[0]);
2851     int j;
2852     for (j=0; j<num_of_clusters; j++)
2853       if (s[j].lab==lab || s[j].max==0) break;
2854     s[j].lab=lab;
2855     s[j].max++;
2856   }
2857   
2858   int max=0;
2859   for (i=0; i<num_of_clusters; i++) 
2860     if (s[i].max>max) {max=s[i].max; lab=s[i].lab;}
2861     
2862   delete[] s;
2863   
2864   for (i=0; i<num_of_clusters; i++) {
2865     AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(i);
2866     if (TMath::Abs(c->fTracks[1]) == lab ||
2867         TMath::Abs(c->fTracks[2]) == lab ) max++;
2868   }
2869   
2870   if (1.-float(max)/num_of_clusters > 0.10) return -lab;
2871   
2872   int tail=int(0.08*nrows);
2873   if (num_of_clusters < tail) return lab;
2874   
2875   max=0;
2876   for (i=1; i<=tail; i++) {
2877     AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(num_of_clusters-i);
2878     if (lab == TMath::Abs(c->fTracks[0]) ||
2879         lab == TMath::Abs(c->fTracks[1]) ||
2880         lab == TMath::Abs(c->fTracks[2])) max++;
2881   }
2882   if (max < int(0.5*tail)) return -lab;
2883   
2884   return lab;
2885 }
2886
2887 //_____________________________________________________________________________
2888 void AliTPCtrack::GetPxPyPz(Double_t& px, Double_t& py, Double_t& pz) const 
2889 {
2890   //-----------------------------------------------------------------
2891   // This function returns reconstructed track momentum in the global system.
2892   //
2893   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2894   //-----------------------------------------------------------------
2895   Double_t pt=TMath::Abs(GetPt()); // GeV/c
2896   Double_t r=x(2)*fX-x(3);
2897   Double_t y0=x(0) + sqrt(1.- r*r)/x(2);
2898   px=-pt*(x(0)-y0)*x(2);    //cos(phi);
2899   py=-pt*(x(3)-fX*x(2));   //sin(phi);
2900   pz=pt*x(4);
2901   Double_t tmp=px*TMath::Cos(fAlpha) - py*TMath::Sin(fAlpha);
2902   py=px*TMath::Sin(fAlpha) + py*TMath::Cos(fAlpha);
2903   px=tmp;  
2904 }
2905
2906 //_____________________________________________________________________________
2907 Double_t AliTPCtrack::GetdEdX(Double_t low, Double_t up) const {
2908   //-----------------------------------------------------------------
2909   // This funtion calculates dE/dX within the "low" and "up" cuts.
2910   //
2911   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2912   //-----------------------------------------------------------------
2913   int ncl=fClusters.GetEntriesFast();
2914   int n=0;
2915   Double_t *q=new Double_t[ncl];
2916   int i;
2917   for (i=0; i<ncl; i++) {
2918      AliTPCcluster *cl=(AliTPCcluster*)(fClusters.UncheckedAt(i));
2919      //     if (cl->fdEdX > 3000) continue;
2920      if (cl->fdEdX <= 0) continue;
2921      q[n++]=cl->fdEdX;
2922   }
2923
2924   //stupid sorting
2925   int swap;
2926   do {
2927     swap=0;
2928     for (i=0; i<n-1; i++) {
2929       if (q[i]<=q[i+1]) continue;
2930       Double_t tmp=q[i]; q[i]=q[i+1]; q[i+1]=tmp;
2931       swap++;
2932     }
2933   } while (swap);
2934
2935   int nl=int(low*n), nu=int(up *n);
2936   Double_t dedx=0.;
2937   for (i=nl; i<=nu; i++) dedx += q[i];
2938   dedx /= (nu-nl+1);
2939   return dedx;
2940 }
2941
2942 //_________________________________________________________________________
2943 //
2944 // Classes for internal tracking use
2945 //_________________________________________________________________________
2946 void AliTPCRow::InsertCluster(const AliTPCcluster* c) {
2947   //-----------------------------------------------------------------------
2948   // Insert a cluster into this pad row in accordence with its y-coordinate
2949   //
2950   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2951   //-----------------------------------------------------------------------
2952   if (num_of_clusters==MAX_CLUSTER_PER_ROW) {
2953     cerr<<"AliTPCRow::InsertCluster(): Too many clusters !\n"; return;
2954   }
2955   if (num_of_clusters==0) {clusters[num_of_clusters++]=c; return;}
2956   int i=Find(c->fY);
2957   memmove(clusters+i+1 ,clusters+i,(num_of_clusters-i)*sizeof(AliTPCcluster*));
2958   clusters[i]=c; num_of_clusters++;
2959 }
2960
2961 int AliTPCRow::Find(Double_t y) const {
2962   //-----------------------------------------------------------------------
2963   // Return the index of the nearest cluster 
2964   //
2965   // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2966   //-----------------------------------------------------------------------
2967   if (y <= clusters[0]->fY) return 0;
2968   if (y > clusters[num_of_clusters-1]->fY) return num_of_clusters;
2969   int b=0, e=num_of_clusters-1, m=(b+e)/2;
2970   for (; b<e; m=(b+e)/2) {
2971     if (y > clusters[m]->fY) b=m+1;
2972     else e=m; 
2973   }
2974   return m;
2975 }
2976