1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
20 ///////////////////////////////////////////////////////////////////////////////
22 // Time Projection Chamber //
23 // This class contains the basic functions for the Time Projection Chamber //
24 // detector. Functions specific to one particular geometry are //
25 // contained in the derived classes //
29 <img src="picts/AliTPCClass.gif">
34 ///////////////////////////////////////////////////////////////////////////////
40 #include <TGeometry.h>
43 #include <TObjectTable.h>
44 #include "TParticle.h"
52 #include "AliTPCParam.h"
54 #include "AliTPCPRF2D.h"
55 #include "AliTPCRF1D.h"
61 //_____________________________________________________________________________
65 // Default constructor
76 fDigParam= new AliTPCD();
77 fDigits = fDigParam->GetArray();
80 //_____________________________________________________________________________
81 AliTPC::AliTPC(const char *name, const char *title)
82 : AliDetector(name,title)
85 // Standard constructor
89 // Initialise arrays of hits and digits
90 fHits = new TClonesArray("AliTPChit", 176);
91 // fDigits = new TClonesArray("AliTPCdigit",10000);
93 fDigParam= new AliTPCD;
94 fDigits = fDigParam->GetArray();
96 AliTPCParam *fTPCParam = &(fDigParam->GetParam());
99 // Initialise counters
103 fNsectors = fTPCParam->GetNSector();
106 fDigitsIndex = new Int_t[fNsectors+1];
107 fClustersIndex = new Int_t[fNsectors+1];
111 // Initialise color attributes
112 SetMarkerColor(kYellow);
115 //_____________________________________________________________________________
127 if (fDigitsIndex) delete [] fDigitsIndex;
128 if (fClustersIndex) delete [] fClustersIndex;
131 //_____________________________________________________________________________
132 void AliTPC::AddCluster(Float_t *hits, Int_t *tracks)
135 // Add a simulated cluster to the list
137 if(!fClusters) fClusters=new TClonesArray("AliTPCcluster",10000);
138 TClonesArray &lclusters = *fClusters;
139 new(lclusters[fNclusters++]) AliTPCcluster(hits,tracks);
142 //_____________________________________________________________________________
143 void AliTPC::AddCluster(const AliTPCcluster &c)
146 // Add a simulated cluster copy to the list
148 if(!fClusters) fClusters=new TClonesArray("AliTPCcluster",10000);
149 TClonesArray &lclusters = *fClusters;
150 new(lclusters[fNclusters++]) AliTPCcluster(c);
153 //_____________________________________________________________________________
154 void AliTPC::AddDigit(Int_t *tracks, Int_t *digits)
157 // Add a TPC digit to the list
159 // TClonesArray &ldigits = *fDigits;
161 TClonesArray &ldigits = *fDigParam->GetArray();
162 new(ldigits[fNdigits++]) AliTPCdigit(tracks,digits);
165 //_____________________________________________________________________________
166 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
169 // Add a hit to the list
171 TClonesArray &lhits = *fHits;
172 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
175 //_____________________________________________________________________________
176 void AliTPC::AddTrack(Float_t *hits)
179 // Add a track to the list of tracks
181 TClonesArray <racks = *fTracks;
182 new(ltracks[fNtracks++]) AliTPCtrack(hits);
185 //_____________________________________________________________________________
186 void AliTPC::AddTrack(const AliTPCtrack& t)
189 // Add a track copy to the list of tracks
191 if(!fTracks) fTracks=new TClonesArray("AliTPCtrack",10000);
192 TClonesArray <racks = *fTracks;
193 new(ltracks[fNtracks++]) AliTPCtrack(t);
196 //_____________________________________________________________________________
197 void AliTPC::BuildGeometry()
200 // Build TPC ROOT TNode geometry for the event display
205 const int kColorTPC=19;
206 char name[5], title[25];
207 const Double_t kDegrad=TMath::Pi()/180;
208 const Double_t kRaddeg=180./TMath::Pi();
210 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
212 Float_t InnerOpenAngle = fTPCParam->GetInnerAngle();
213 Float_t OuterOpenAngle = fTPCParam->GetOuterAngle();
215 Float_t InnerAngleShift = fTPCParam->GetInnerAngleShift();
216 Float_t OuterAngleShift = fTPCParam->GetOuterAngleShift();
218 Int_t nLo = fTPCParam->GetNInnerSector()/2;
219 Int_t nHi = fTPCParam->GetNOuterSector()/2;
221 const Double_t loAng = (Double_t)TMath::Nint(InnerOpenAngle*kRaddeg);
222 const Double_t hiAng = (Double_t)TMath::Nint(OuterOpenAngle*kRaddeg);
223 const Double_t loAngSh = (Double_t)TMath::Nint(InnerAngleShift*kRaddeg);
224 const Double_t hiAngSh = (Double_t)TMath::Nint(OuterAngleShift*kRaddeg);
227 const Double_t loCorr = 1/TMath::Cos(0.5*loAng*kDegrad);
228 const Double_t hiCorr = 1/TMath::Cos(0.5*hiAng*kDegrad);
234 // Get ALICE top node
237 Top=gAlice->GetGeometry()->GetNode("alice");
241 rl = fTPCParam->GetInSecLowEdge();
242 ru = fTPCParam->GetInSecUpEdge();
246 sprintf(name,"LS%2.2d",i);
248 sprintf(title,"TPC low sector %3d",i);
251 tubs = new TTUBS(name,title,"void",rl*loCorr,ru*loCorr,250.,
252 loAng*(i-0.5)+loAngSh,loAng*(i+0.5)+loAngSh);
253 tubs->SetNumberOfDivisions(1);
255 Node = new TNode(name,title,name,0,0,0,"");
256 Node->SetLineColor(kColorTPC);
262 rl = fTPCParam->GetOuSecLowEdge();
263 ru = fTPCParam->GetOuSecUpEdge();
266 sprintf(name,"US%2.2d",i);
268 sprintf(title,"TPC upper sector %d",i);
270 tubs = new TTUBS(name,title,"void",rl*hiCorr,ru*hiCorr,250,
271 hiAng*(i-0.5)+hiAngSh,hiAng*(i+0.5)+hiAngSh);
272 tubs->SetNumberOfDivisions(1);
274 Node = new TNode(name,title,name,0,0,0,"");
275 Node->SetLineColor(kColorTPC);
282 //_____________________________________________________________________________
283 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
286 // Calculate distance from TPC to mouse on the display
292 //_____________________________________________________________________________
293 static Double_t SigmaY2(Double_t r, Double_t tgl, Double_t pt)
296 // Parametrised error of the cluster reconstruction (pad direction)
298 pt=TMath::Abs(pt)*1000.;
301 Double_t s=a_rphi - b_rphi*r*tgl + c_rphi*x*x + d_rphi*x;
302 if (s<0.4e-3) s=0.4e-3;
303 s*=1.3; //Iouri Belikov
307 //_____________________________________________________________________________
308 static Double_t SigmaZ2(Double_t r, Double_t tgl)
311 // Parametrised error of the cluster reconstruction (drift direction)
314 Double_t s=a_z - b_z*r*tgl + c_z*tgl*tgl;
315 if (s<0.4e-3) s=0.4e-3;
316 s*=1.3; //Iouri Belikov
320 //_____________________________________________________________________________
321 inline Double_t f1(Double_t x1,Double_t y1,
322 Double_t x2,Double_t y2,
323 Double_t x3,Double_t y3)
325 //-----------------------------------------------------------------
326 // Initial approximation of the track curvature
328 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
329 //-----------------------------------------------------------------
330 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
331 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
332 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
333 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
334 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
336 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
338 return -xr*yr/sqrt(xr*xr+yr*yr);
342 //_____________________________________________________________________________
343 inline Double_t f2(Double_t x1,Double_t y1,
344 Double_t x2,Double_t y2,
345 Double_t x3,Double_t y3)
347 //-----------------------------------------------------------------
348 // Initial approximation of the track curvature times center of curvature
350 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
351 //-----------------------------------------------------------------
352 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
353 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
354 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
355 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
356 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
358 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
360 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
363 //_____________________________________________________________________________
364 inline Double_t f3(Double_t x1,Double_t y1,
365 Double_t x2,Double_t y2,
366 Double_t z1,Double_t z2)
368 //-----------------------------------------------------------------
369 // Initial approximation of the tangent of the track dip angle
371 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
372 //-----------------------------------------------------------------
373 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
376 //_____________________________________________________________________________
377 static int FindProlongation(AliTPCtrack& t, const AliTPCSector *sec,
380 //-----------------------------------------------------------------
381 // This function tries to find a track prolongation.
383 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
384 //-----------------------------------------------------------------
385 const int ROWS_TO_SKIP=int(0.5*sec->GetNRows());
386 const Float_t MAX_CHI2=12.;
387 int try_again=ROWS_TO_SKIP;
388 Double_t alpha=sec->GetAlpha();
389 int ns=int(2*TMath::Pi()/alpha+0.5);
391 for (int nr=sec->GetRowNumber(t.GetX())-1; nr>=rf; nr--) {
392 Double_t x=sec->GetX(nr), ymax=sec->GetMaxY(nr);
393 if (!t.PropagateTo(x)) return 0;
396 Double_t max_chi2=MAX_CHI2;
397 const AliTPCRow& row=sec[s][nr];
398 Double_t sy2=SigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
399 Double_t sz2=SigmaZ2(t.GetX(),t.GetTgl());
400 Double_t road=5.*sqrt(t.GetSigmaY2() + sy2), y=t.GetY(), z=t.GetZ();
403 if (t>4) cerr<<t<<" FindProlongation warning: Too broad road !\n";
408 for (int i=row.Find(y-road); i<row; i++) {
409 AliTPCcluster* c=(AliTPCcluster*)(row[i]);
410 if (c->fY > y+road) break;
411 if (c->IsUsed()) continue;
412 if ((c->fZ - z)*(c->fZ - z) > 25.*(t.GetSigmaZ2() + sz2)) continue;
413 Double_t chi2=t.GetPredictedChi2(c);
414 if (chi2 > max_chi2) continue;
420 t.Update(cl,max_chi2);
421 Double_t ll=TMath::Sqrt((1+t.GetTgl()*t.GetTgl())/
422 (1-(t.GetC()*x-t.GetEta())*(t.GetC()*x-t.GetEta())));
423 cl->fdEdX = cl->fQ/ll;
424 try_again=ROWS_TO_SKIP;
426 if (try_again==0) break;
429 if (!t.Rotate(alpha)) return 0;
430 } else if (y <-ymax) {
432 if (!t.Rotate(-alpha)) return 0;
442 //_____________________________________________________________________________
443 static void MakeSeeds(TObjArray& seeds,const AliTPCSector *sec, int max_sec,
446 //-----------------------------------------------------------------
447 // This function creates track seeds.
449 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
450 //-----------------------------------------------------------------
451 TMatrix C(5,5); TVector x(5);
452 double alpha=sec->GetAlpha(), shift=sec->GetAlphaShift();
453 double cs=cos(alpha), sn=sin(alpha);
454 for (int ns=0; ns<max_sec; ns++) {
455 int nl=sec[(ns-1+max_sec)%max_sec][i2];
457 int nu=sec[(ns+1)%max_sec][i2];
458 const AliTPCRow& r1=sec[ns][i1];
459 for (int is=0; is < r1; is++) {
460 double x1=sec->GetX(i1), y1=r1[is]->fY, z1=r1[is]->fZ;
461 for (int js=0; js < nl+nm+nu; js++) {
462 const AliTPCcluster *cl;
464 double x2=sec->GetX(i2), y2, z2, tmp;
467 ks=(ns-1+max_sec)%max_sec;
468 const AliTPCRow& r2=sec[(ns-1+max_sec)%max_sec][i2];
470 y2=cl->fY; z2=cl->fZ;
472 y2 =-x2*sn+y2*cs; x2=tmp;
476 const AliTPCRow& r2=sec[ns][i2];
478 y2=cl->fY; z2=cl->fZ;
481 const AliTPCRow& r2=sec[(ns+1)%max_sec][i2];
483 y2=cl->fY; z2=cl->fZ;
485 y2 =x2*sn+y2*cs; x2=tmp;
488 double d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
489 if (d==0.) {cerr<<"MakeSeeds warning: Straight seed !\n"; continue;}
493 x(2)=f1(x1,y1,x2,y2,0.,0.);
494 x(3)=f2(x1,y1,x2,y2,0.,0.);
495 x(4)=f3(x1,y1,x2,y2,z1,z2);
497 if (TMath::Abs(x(2)*x1-x(3)) >= 0.999) continue;
499 if (TMath::Abs(x(4)) > 1.2) continue;
501 Double_t a=asin(x(3));
502 Double_t zv=z1 - x(4)/x(2)*(a+asin(x(2)*x1-x(3)));
503 if (TMath::Abs(zv)>33.) continue;
505 TMatrix X(6,6); X=0.;
506 X(0,0)=r1[is]->fSigmaY2; X(1,1)=r1[is]->fSigmaZ2;
507 X(2,2)=cl->fSigmaY2; X(3,3)=cl->fSigmaZ2;
508 X(4,4)=3./12.; X(5,5)=3./12.;
509 TMatrix F(5,6); F.UnitMatrix();
510 Double_t sy=sqrt(X(0,0)), sz=sqrt(X(1,1));
511 F(2,0)=(f1(x1,y1+sy,x2,y2,0.,0.)-x(2))/sy;
512 F(2,2)=(f1(x1,y1,x2,y2+sy,0.,0.)-x(2))/sy;
513 F(2,4)=(f1(x1,y1,x2,y2,0.,0.+sy)-x(2))/sy;
514 F(3,0)=(f2(x1,y1+sy,x2,y2,0.,0.)-x(3))/sy;
515 F(3,2)=(f2(x1,y1,x2,y2+sy,0.,0.)-x(3))/sy;
516 F(3,4)=(f2(x1,y1,x2,y2,0.,0.+sy)-x(3))/sy;
517 F(4,0)=(f3(x1,y1+sy,x2,y2,z1,z2)-x(4))/sy;
518 F(4,1)=(f3(x1,y1,x2,y2,z1+sz,z2)-x(4))/sz;
519 F(4,2)=(f3(x1,y1,x2,y2+sy,z1,z2)-x(4))/sy;
520 F(4,3)=(f3(x1,y1,x2,y2,z1,z2+sz)-x(4))/sz;
524 TMatrix t(F,TMatrix::kMult,X);
525 C.Mult(t,TMatrix(TMatrix::kTransposed,F));
527 AliTPCtrack *track=new AliTPCtrack(r1[is], x, C, x1, ns*alpha+shift);
528 int rc=FindProlongation(*track,sec,ns,i2);
529 if (rc<0 || *track<(i1-i2)/2) delete track;
530 else seeds.AddLast(track);
536 //_____________________________________________________________________________
537 AliTPCParam *AliTPCSector::param;
538 void AliTPC::Clusters2Tracks()
540 //-----------------------------------------------------------------
541 // This is a track finder.
543 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
544 //-----------------------------------------------------------------
545 if (!fClusters) return;
547 AliTPCParam *p=&fDigParam->GetParam();
548 AliTPCSector::SetParam(p);
550 const int nis=p->GetNInnerSector()/2;
551 AliTPCSSector *ssec=new AliTPCSSector[nis];
552 int nrow_low=ssec->GetNRows();
554 const int nos=p->GetNOuterSector()/2;
555 AliTPCLSector *lsec=new AliTPCLSector[nos];
556 int nrow_up=lsec->GetNRows();
558 int ncl=fClusters->GetEntriesFast();
560 AliTPCcluster *c=(AliTPCcluster*)fClusters->UncheckedAt(ncl);
561 Int_t sec=c->fSector, row=c->fPadRow;
563 ssec[sec%nis][row].InsertCluster(c);
566 lsec[sec%nos][row].InsertCluster(c);
570 TObjArray seeds(20000);
572 int nrows=nrow_low+nrow_up;
573 int gap=int(0.125*nrows), shift=int(0.5*gap);
574 MakeSeeds(seeds, lsec, nos, nrow_up-1, nrow_up-1-gap);
575 MakeSeeds(seeds, lsec, nos, nrow_up-1-shift, nrow_up-1-shift-gap);
580 int nseed=seeds.GetEntriesFast();
582 for (int s=0; s<nseed; s++) {
583 AliTPCtrack& t=*((AliTPCtrack*)seeds.UncheckedAt(s));
584 double alpha=t.GetAlpha();
585 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
586 if (alpha < 0. ) alpha += 2.*TMath::Pi();
587 int ns=int(alpha/lsec->GetAlpha())%nos;
589 if (!FindProlongation(t,lsec,ns)) continue;
591 alpha=t.GetAlpha() + 0.5*ssec->GetAlpha() - ssec->GetAlphaShift();
592 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
593 if (alpha < 0. ) alpha += 2.*TMath::Pi();
594 ns=int(alpha/ssec->GetAlpha())%nis; //index of the inner sector needed
596 alpha=ns*ssec->GetAlpha() - t.GetAlpha();
597 if (!t.Rotate(alpha)) continue;
599 if (!FindProlongation(t,ssec,ns)) continue;
601 if (t < int(0.4*nrows)) continue;
612 //_____________________________________________________________________________
613 void AliTPC::CreateMaterials()
615 //-----------------------------------------------
616 // Create Materials for for TPC
617 //-----------------------------------------------
619 //-----------------------------------------------------------------
620 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
621 //-----------------------------------------------------------------
623 Int_t ISXFLD=gAlice->Field()->Integ();
624 Float_t SXMGMX=gAlice->Field()->Max();
626 Float_t amat[5]; // atomic numbers
627 Float_t zmat[5]; // z
628 Float_t wmat[5]; // proportions
632 // ********************* Gases *******************
634 //--------------------------------------------------------------
636 //--------------------------------------------------------------
641 Float_t a_ne = 20.18;
646 AliMaterial(20,"Ne",a_ne,z_ne,density,999.,999.);
650 Float_t a_ar = 39.948;
655 AliMaterial(21,"Ar",a_ar,z_ar,density,999.,999.);
663 //--------------------------------------------------------------
665 //--------------------------------------------------------------
682 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
684 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
699 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
701 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
716 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
718 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
720 //----------------------------------------------------------------
721 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
722 //----------------------------------------------------------------
730 Float_t a,z,rho,absl,X0,buf[1];
733 for(nc = 0;nc<fNoComp;nc++)
736 // retrive material constants
738 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
743 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
745 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? a_pure[nnc] : amol[nnc]);
746 density += fMixtProp[nc]*rho; // density of the mixture
750 // mixture proportions by weight!
752 for(nc = 0;nc<fNoComp;nc++)
755 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
757 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ? a_pure[nnc] : amol[nnc])/am;
761 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
762 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
763 AliMixture(33,"Drift gas 3",amat,zmat,density,fNoComp,wmat);
765 AliMedium(2, "Drift gas 1", 31, 0, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
766 AliMedium(3, "Drift gas 2", 32, 0, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
767 AliMedium(4, "Drift gas 3", 33, 1, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
771 AliMaterial(24, "Air", 14.61, 7.3, .001205, 30420., 67500.);
773 AliMedium(24, "Air", 24, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
775 //----------------------------------------------------------------------
777 //----------------------------------------------------------------------
781 AliMaterial(30, "Al", 26.98, 13., 2.7, 8.9, 37.2);
783 AliMedium(0, "Al",30, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
787 AliMaterial(31, "Si", 28.086, 14.,2.33, 9.36, 999.);
789 AliMedium(7, "Al",31, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
808 AliMixture(32, "Mylar",amat,zmat,density,-3,wmat);
810 AliMedium(5, "Mylar",32, 0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
817 AliMaterial(33,"C normal",12.011,6.,2.265,18.8,999.);
819 AliMedium(6,"C normal",33,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
821 // G10 for inner and outr field cage
822 // G10 is 60% SiO2 + 40% epoxy, right now I use A and Z for SiO2
838 AliMixture(34,"G10 aux.",amat,zmat,density,-2,wmat);
841 gMC->Gfmate((*fIdmate)[34],namate,a,z,rho,X0,absl,buf,nbuf);
843 Float_t thickX0 = 0.0052; // field cage in X0 units
845 Float_t thick = 2.; // in cm
849 rhoFactor = X0*thickX0/thick;
850 density = rho*rhoFactor;
852 AliMaterial(35,"G10-fc",a,z,density,999.,999.);
854 AliMedium(8,"G10-fc",35,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
856 thickX0 = 0.0027; // inner vessel (eta <0.9)
858 rhoFactor = X0*thickX0/thick;
859 density = rho*rhoFactor;
861 AliMaterial(36,"G10-iv",a,z,density,999.,999.);
863 AliMedium(9,"G10-iv",36,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
867 gMC->Gfmate((*fIdmate)[33],namate,a,z,rho,X0,absl,buf,nbuf);
869 thickX0 = 0.0133; // outer vessel
871 rhoFactor = X0*thickX0/thick;
872 density = rho*rhoFactor;
875 AliMaterial(37,"C-ov",a,z,density,999.,999.);
877 AliMedium(10,"C-ov",37,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
879 thickX0=0.015; // inner vessel (cone, eta > 0.9)
881 rhoFactor = X0*thickX0/thick;
882 density = rho*rhoFactor;
884 AliMaterial(38,"C-ivc",a,z,density,999.,999.);
886 AliMedium(11,"C-ivc",38,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
890 AliMedium(12,"CO2",10,0, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
896 //_____________________________________________________________________________
898 const AliTPCdigit *dig;
900 Bin() {dig=0; idx=-1;}
903 struct PreCluster : public AliTPCcluster {
904 const AliTPCdigit* summit; //pointer to the maximum digit of this precluster
905 int idx; //index in AliTPC::fClusters
906 int npeaks; //number of peaks in this precluster
907 int ndigits; //number of digits in this precluster
910 PreCluster::PreCluster() : AliTPCcluster() {npeaks=ndigits=0;}
912 //_____________________________________________________________________________
913 static void FindPreCluster(int i,int j,int maxj,Bin *bins,PreCluster &c)
915 //-----------------------------------------------------------------
916 // This function looks for "preclusters".
918 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
919 //-----------------------------------------------------------------
920 Bin& b=bins[i*maxj+j];
921 double q=(double)TMath::Abs(b.dig->fSignal);
923 if (b.idx >= 0 && b.idx != c.idx) {
928 if (q > TMath::Abs(c.summit->fSignal)) c.summit=b.dig;
937 b.dig = 0; b.idx = c.idx;
939 if (bins[(i-1)*maxj+j].dig) FindPreCluster(i-1,j,maxj,bins,c);
940 if (bins[i*maxj+(j-1)].dig) FindPreCluster(i,j-1,maxj,bins,c);
941 if (bins[(i+1)*maxj+j].dig) FindPreCluster(i+1,j,maxj,bins,c);
942 if (bins[i*maxj+(j+1)].dig) FindPreCluster(i,j+1,maxj,bins,c);
946 //_____________________________________________________________________________
947 void AliTPC::Digits2Clusters()
949 //-----------------------------------------------------------------
950 // This is a simple cluster finder.
952 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
953 //-----------------------------------------------------------------
954 AliTPCParam *par = &(fDigParam->GetParam());
956 int inp=par->GetNPads(0, par->GetNRowLow()-1);
957 int onp=par->GetNPads(par->GetNSector()-1,par->GetNRowUp() -1);
958 const int MAXY=(inp>onp) ? inp+2 : onp+2;
959 const int MAXTBKT=int((z_end+6*par->GetZSigma())/par->GetZWidth())+1;
960 const int MAXZ=MAXTBKT+2;
961 const int THRESHOLD=20;
963 TTree *t=(TTree*)gDirectory->Get("TreeD0_Param1");
964 t->GetBranch("Digits")->SetAddress(&fDigits);
965 Int_t sectors_by_rows=(Int_t)t->GetEntries();
969 for (Int_t n=0; n<sectors_by_rows; n++) {
970 if (!t->GetEvent(n)) continue;
971 Bin *bins=new Bin[MAXY*MAXZ];
972 AliTPCdigit *dig=(AliTPCdigit*)fDigits->UncheckedAt(0);
973 int sec=dig->fSector, row=dig->fPadRow;
974 int ndigits=fDigits->GetEntriesFast();
978 int nis=par->GetNInnerSector(), nos=par->GetNOuterSector();
980 npads = par->GetNPadsLow(row);
981 sign = (sec < nis/2) ? 1 : -1;
983 npads = par->GetNPadsUp(row);
984 sign = ((sec-nis) < nos/2) ? 1 : -1;
989 for (ndig=0; ndig<ndigits; ndig++) {
990 dig=(AliTPCdigit*)fDigits->UncheckedAt(ndig);
991 int i=dig->fPad+1, j=dig->fTime+1;
993 cerr<<"AliTPC::Digits2Clusters error: pad number is out of range ! ";
994 cerr<<i<<' '<<npads<<endl;
998 cerr<<"AliTPC::Digits2Clusters error: time bucket is out of range ! ";
999 cerr<<j<<' '<<MAXTBKT<<endl;
1002 if (dig->fSignal >= THRESHOLD) bins[i*MAXZ+j].dig=dig;
1003 if (i==1 || i==npads || j==1 || j==MAXTBKT) dig->fSignal*=-1;
1009 for (i=1; i<MAXY-1; i++) {
1010 for (j=1; j<MAXZ-1; j++) {
1011 if (bins[i*MAXZ+j].dig == 0) continue;
1012 PreCluster c; c.summit=bins[i*MAXZ+j].dig; c.idx=ncls;
1013 FindPreCluster(i, j, MAXZ, bins, c);
1017 double s2 = c.fSigmaY2/c.fQ - c.fY*c.fY;
1018 c.fSigmaY2 = s2 + 1./12.;
1019 c.fSigmaY2 *= par->GetPadPitchWidth()*par->GetPadPitchWidth();
1020 if (s2 != 0.) c.fSigmaY2 *= 0.17;
1022 s2 = c.fSigmaZ2/c.fQ - c.fZ*c.fZ;
1023 c.fSigmaZ2 = s2 + 1./12.;
1024 c.fSigmaZ2 *= par->GetZWidth()*par->GetZWidth();
1025 if (s2 != 0.) c.fSigmaZ2 *= 0.41;
1027 c.fY = (c.fY - 0.5 - 0.5*npads)*par->GetPadPitchWidth();
1028 c.fZ = par->GetZWidth()*c.fZ;
1029 c.fZ -= 3.*par->GetZSigma(); // PASA delay
1030 c.fZ = sign*(z_end - c.fZ);
1034 c.fTracks[0]=c.summit->fTracks[0];
1035 c.fTracks[1]=c.summit->fTracks[1];
1036 c.fTracks[2]=c.summit->fTracks[2];
1038 if (c.summit->fSignal<0) {
1043 AddCluster(c); ncls++; ncl++;
1047 for (ndig=0; ndig<ndigits; ndig++) {
1048 dig=(AliTPCdigit*)fDigits->UncheckedAt(ndig);
1049 int i=dig->fPad+1, j=dig->fTime+1;
1051 cerr<<"AliTPC::Digits2Clusters error: pad number is out of range ! ";
1052 cerr<<i<<' '<<npads<<endl;
1056 cerr<<"AliTPC::Digits2Clusters error: time bucket is out of range ! ";
1057 cerr<<j<<' '<<MAXTBKT<<endl;
1060 if (TMath::Abs(dig->fSignal)>=par->GetZeroSup()) bins[i*MAXZ+j].dig=dig;
1063 for (i=1; i<MAXY-1; i++) {
1064 for (j=1; j<MAXZ-1; j++) {
1065 if (bins[i*MAXZ+j].dig == 0) continue;
1066 PreCluster c; c.summit=bins[i*MAXZ+j].dig; c.idx=ncls;
1067 FindPreCluster(i, j, MAXZ, bins, c);
1068 if (c.ndigits < 2) continue; //noise cluster
1069 if (c.npeaks>1) continue; //overlapped cluster
1073 double s2 = c.fSigmaY2/c.fQ - c.fY*c.fY;
1074 c.fSigmaY2 = s2 + 1./12.;
1075 c.fSigmaY2 *= par->GetPadPitchWidth()*par->GetPadPitchWidth();
1076 if (s2 != 0.) c.fSigmaY2 *= 0.04;
1078 s2 = c.fSigmaZ2/c.fQ - c.fZ*c.fZ;
1079 c.fSigmaZ2 = s2 + 1./12.;
1080 c.fSigmaZ2 *= par->GetZWidth()*par->GetZWidth();
1081 if (s2 != 0.) c.fSigmaZ2 *= 0.10;
1083 c.fY = (c.fY - 0.5 - 0.5*npads)*par->GetPadPitchWidth();
1084 c.fZ = par->GetZWidth()*c.fZ;
1085 c.fZ -= 3.*par->GetZSigma(); // PASA delay
1086 c.fZ = sign*(z_end - c.fZ);
1090 c.fTracks[0]=c.summit->fTracks[0];
1091 c.fTracks[1]=c.summit->fTracks[1];
1092 c.fTracks[2]=c.summit->fTracks[2];
1094 if (c.summit->fSignal<0) {
1099 if (c.npeaks==0) {AddCluster(c); ncls++; ncl++;}
1101 new ((*fClusters)[c.idx]) AliTPCcluster(c);
1106 cerr<<"sector, row, digits, clusters: "
1107 <<sec<<' '<<row<<' '<<ndigits<<' '<<ncl<<" \r";
1115 //_____________________________________________________________________________
1116 void AliTPC::ElDiff(Float_t *xyz)
1118 //--------------------------------------------------
1119 // calculates the diffusion of a single electron
1120 //--------------------------------------------------
1122 //-----------------------------------------------------------------
1123 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1124 //-----------------------------------------------------------------
1125 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1130 driftl=z_end-TMath::Abs(xyz[2]);
1132 if(driftl<0.01) driftl=0.01;
1134 // check the attachment
1136 driftl=TMath::Sqrt(driftl);
1138 // Float_t sig_t = driftl*diff_t;
1139 //Float_t sig_l = driftl*diff_l;
1140 Float_t sig_t = driftl*fTPCParam->GetDiffT();
1141 Float_t sig_l = driftl*fTPCParam->GetDiffL();
1142 xyz[0]=gRandom->Gaus(xyz[0],sig_t);
1143 xyz[1]=gRandom->Gaus(xyz[1],sig_t);
1144 xyz[2]=gRandom->Gaus(xyz[2],sig_l);
1146 if (TMath::Abs(xyz[2])>z_end){
1147 xyz[2]=TMath::Sign(z_end,z0);
1150 Float_t eps = 0.0001;
1151 xyz[2]=TMath::Sign(eps,z0);
1155 //_____________________________________________________________________________
1156 void AliTPC::Hits2Clusters()
1158 //--------------------------------------------------------
1159 // TPC simple cluster generator from hits
1160 // obtained from the TPC Fast Simulator
1161 // The point errors are taken from the parametrization
1162 //--------------------------------------------------------
1164 //-----------------------------------------------------------------
1165 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1166 //-----------------------------------------------------------------
1168 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1169 Float_t sigma_rphi,sigma_z,cl_rphi,cl_z;
1171 TParticle *particle; // pointer to a given particle
1172 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1173 TClonesArray *Particles; //pointer to the particle list
1177 Float_t pl,pt,tanth,rpad,ratio;
1180 //---------------------------------------------------------------
1181 // Get the access to the tracks
1182 //---------------------------------------------------------------
1184 TTree *TH = gAlice->TreeH();
1185 Stat_t ntracks = TH->GetEntries();
1186 Particles=gAlice->Particles();
1188 //------------------------------------------------------------
1189 // Loop over all sectors (72 sectors)
1190 // Sectors 0-35 are lower sectors, 0-17 z>0, 18-35 z<0
1191 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1193 // First cluster for sector 0 starts at "0"
1194 //------------------------------------------------------------
1197 //fClustersIndex[0] = 0;
1200 int Nsectors=fDigParam->GetParam().GetNSector();
1201 for(Int_t isec=0; isec<Nsectors; isec++){
1203 fTPCParam->AdjustAngles(isec,cph,sph);
1205 //------------------------------------------------------------
1207 //------------------------------------------------------------
1209 for(Int_t track=0;track<ntracks;track++){
1211 TH->GetEvent(track);
1213 // Get number of the TPC hits and a pointer
1216 nhits=fHits->GetEntriesFast();
1220 for(Int_t hit=0;hit<nhits;hit++){
1221 tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
1222 sector=tpcHit->fSector; // sector number
1223 if(sector != isec) continue; //terminate iteration
1224 ipart=tpcHit->fTrack;
1225 particle=(TParticle*)Particles->UncheckedAt(ipart);
1228 if(pt < 1.e-9) pt=1.e-9;
1230 tanth = TMath::Abs(tanth);
1231 rpad=TMath::Sqrt(tpcHit->fX*tpcHit->fX + tpcHit->fY*tpcHit->fY);
1232 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
1234 // space-point resolutions
1236 sigma_rphi=SigmaY2(rpad,tanth,pt);
1237 sigma_z =SigmaZ2(rpad,tanth );
1241 cl_rphi=ac_rphi-bc_rphi*rpad*tanth+cc_rphi*ratio*ratio;
1242 cl_z=ac_z-bc_z*rpad*tanth+cc_z*tanth*tanth;
1244 // temporary protection
1246 if(sigma_rphi < 0.) sigma_rphi=0.4e-3;
1247 if(sigma_z < 0.) sigma_z=0.4e-3;
1248 if(cl_rphi < 0.) cl_rphi=2.5e-3;
1249 if(cl_z < 0.) cl_z=2.5e-5;
1252 // smearing --> rotate sectors firstly,
1253 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1255 Float_t xprim= tpcHit->fX*cph + tpcHit->fY*sph;
1256 Float_t yprim=-tpcHit->fX*sph + tpcHit->fY*cph;
1257 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigma_rphi)); // y
1258 Double_t alpha=(sector < fTPCParam->GetNInnerSector()) ?
1259 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1260 if (TMath::Abs(xyz[0]/xprim) > TMath::Tan(0.5*alpha)) xyz[0]=yprim;
1261 xyz[1]=gRandom->Gaus(tpcHit->fZ,TMath::Sqrt(sigma_z)); // z
1262 if (TMath::Abs(xyz[1]) > 250) xyz[1]=tpcHit->fZ;
1263 xyz[2]=tpcHit->fQ+1;// q; let it be not equal to zero (Y.Belikov)
1264 xyz[3]=sigma_rphi; // fSigmaY2
1265 xyz[4]=sigma_z; // fSigmaZ2
1267 // and finally add the cluster
1268 Int_t tracks[5]={tpcHit->fTrack, -1, -1, sector, tpcHit->fPadRow};
1269 AddCluster(xyz,tracks);
1271 } // end of loop over hits
1272 } // end of loop over tracks
1274 //fClustersIndex[isec] = fNclusters; // update clusters index
1276 } // end of loop over sectors
1278 //fClustersIndex[fNsectors]--; // set end of the clusters buffer
1283 void AliTPC::Hits2Digits()
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
1291 int Nsectors=fDigParam->GetParam().GetNSector();
1292 for(Int_t isec=0;isec<Nsectors;isec++) Hits2DigitsSector(isec);
1296 //_____________________________________________________________________________
1297 void AliTPC::Hits2DigitsSector(Int_t isec)
1299 //-------------------------------------------------------------------
1300 // TPC conversion from hits to digits.
1301 //-------------------------------------------------------------------
1303 //-----------------------------------------------------------------
1304 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1305 //-----------------------------------------------------------------
1307 //-------------------------------------------------------
1308 // Get the access to the track hits
1309 //-------------------------------------------------------
1311 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1312 TTree *TH = gAlice->TreeH(); // pointer to the hits tree
1313 Stat_t ntracks = TH->GetEntries();
1317 //-------------------------------------------
1318 // Only if there are any tracks...
1319 //-------------------------------------------
1322 // TObjArrays for three neighouring pad-rows
1324 TObjArray **rowTriplet = new TObjArray* [3];
1326 // TObjArray-s for each pad-row
1330 for (Int_t trip=0;trip<3;trip++){
1331 rowTriplet[trip]=new TObjArray;
1336 printf("*** Processing sector number %d ***\n",isec);
1338 Int_t nrows =fTPCParam->GetNRow(isec);
1340 row= new TObjArray* [nrows];
1342 MakeSector(isec,nrows,TH,ntracks,row);
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 //--------------------------------------------------------
1353 for (i=0;i<nrows;i++){
1355 // Triplets for i = 0 and i=1 are identical!
1356 // The same for i = nrows-1 and nrows!
1358 if(i != 1 && i != nrows-1){
1359 MakeTriplet(i,rowTriplet,row);
1362 DigitizeRow(i,isec,rowTriplet);
1366 Int_t ndig=fDigParam->GetArray()->GetEntriesFast();
1368 printf("*** Sector, row, digits %d %d %d ***\n",isec,i,ndig);
1370 ResetDigits(); // reset digits for this row after storing them
1372 } // end of the sector digitization
1374 // delete the last triplet
1376 for (i=0;i<3;i++) rowTriplet[i]->Delete();
1378 delete [] row; // delete the array of pointers to TObjArray-s
1381 } // end of Hits2Digits
1382 //_____________________________________________________________________________
1383 void AliTPC::MakeTriplet(Int_t row,
1384 TObjArray **rowTriplet, TObjArray **prow)
1386 //------------------------------------------------------------------
1387 // Makes the "triplet" of the neighbouring pad-row for the
1388 // digitization including the cross-talk between the pad-rows
1389 //------------------------------------------------------------------
1391 //-----------------------------------------------------------------
1392 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1393 //-----------------------------------------------------------------
1395 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1396 Float_t gasgain = fTPCParam->GetGasGain();
1399 Int_t nElements,nElectrons;
1404 //-------------------------------------------------------------------
1405 // pv is an "old" track, i.e. label + triplets of (x,y,z)
1406 // for each electron
1408 //-------------------------------------------------------------------
1414 if(row == 0 || row == 1){
1416 // create entire triplet for the first AND the second row
1418 nTracks[0] = prow[0]->GetEntries();
1419 nTracks[1] = prow[1]->GetEntries();
1420 nTracks[2] = prow[2]->GetEntries();
1422 for(i1=0;i1<3;i1++){
1423 nt = nTracks[i1]; // number of tracks for this row
1425 for(i2=0;i2<nt;i2++){
1426 pv = (TVector*)prow[i1]->At(i2);
1428 nElements = pv->GetNrows();
1429 nElectrons = (nElements-1)/3;
1431 tv = new TVector(4*nElectrons+1); // create TVector for a modified track
1433 v2(0)=v1(0); //track label
1435 for(nel=0;nel<nElectrons;nel++){
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
1446 // Add this track to a row
1449 rowTriplet[i1]->Add(tv);
1452 } // end of loop over tracks for this row
1454 prow[i1]->Delete(); // remove "old" tracks
1455 delete prow[i1]; // delete TObjArray for this row
1456 prow[i1]=0; // set pointer to NULL
1458 } // end of loop over row triplets
1464 rowTriplet[0]->Delete(); // remove old lower row
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
1471 //-------------------------------------------
1472 // shift new tracks downwards
1473 //-------------------------------------------
1475 for(i1=0;i1<nTracks[0];i1++){
1476 pv=(TVector*)rowTriplet[1]->At(i1);
1477 rowTriplet[0]->Add(pv);
1480 rowTriplet[1]->Clear(); // leave tracks on the heap!!!
1482 for(i1=0;i1<nTracks[1];i1++){
1483 pv=(TVector*)rowTriplet[2]->At(i1);
1484 rowTriplet[1]->Add(pv);
1487 rowTriplet[2]->Clear(); // leave tracks on the heap!!!
1489 //---------------------------------------------
1490 // Create new upper row
1491 //---------------------------------------------
1495 for(i1=0;i1<nTracks[2];i1++){
1496 pv = (TVector*)prow[row+1]->At(i1);
1498 nElements = pv->GetNrows();
1499 nElectrons = (nElements-1)/3;
1501 tv = new TVector(4*nElectrons+1); // create TVector for a modified track
1503 v2(0)=v1(0); //track label
1505 for(nel=0;nel<nElectrons;nel++){
1509 // Avalanche, including fluctuations
1510 Int_t aval = (Int_t)
1511 (-gasgain*TMath::Log(gRandom->Rndm()));
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
1519 rowTriplet[2]->Add(tv);
1521 } // end of loop over tracks
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
1529 } // end of MakeTriplet
1530 //_____________________________________________________________________________
1531 void AliTPC::ExB(Float_t *xyz)
1533 //------------------------------------------------
1534 // ExB at the wires and wire number calulation
1535 //------------------------------------------------
1537 //-----------------------------------------------------------------
1538 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1539 //-----------------------------------------------------------------
1540 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1543 fTPCParam->GetWire(x1); //calculate nearest wire position
1544 Float_t dx=xyz[0]-x1;
1545 xyz[1]+=dx*fTPCParam->GetOmegaTau();
1548 //_____________________________________________________________________________
1549 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rowTriplet)
1551 //-----------------------------------------------------------
1552 // Single row digitization, coupling from the neighbouring
1553 // rows taken into account
1554 //-----------------------------------------------------------
1556 //-----------------------------------------------------------------
1557 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1558 //-----------------------------------------------------------------
1561 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1562 Float_t chipgain= fTPCParam->GetChipGain();
1563 Float_t zerosup = fTPCParam->GetZeroSup();
1564 Int_t nrows =fTPCParam->GetNRow(isec);
1566 int((z_end+6*fTPCParam->GetZSigma())/fTPCParam->GetZWidth())+1;
1570 Int_t IndexRange[4];
1575 // iFlag = 0 -> inner row, iFlag = 1 -> the middle one, iFlag = 2 -> the outer one
1578 nTracks[0]=rowTriplet[0]->GetEntries(); // lower row
1579 nTracks[1]=rowTriplet[1]->GetEntries(); // middle row
1580 nTracks[2]=rowTriplet[2]->GetEntries(); // upper row
1585 n_of_pads[0]=fTPCParam->GetNPads(isec,0);
1586 n_of_pads[1]=fTPCParam->GetNPads(isec,1);
1588 else if(irow == nrows-1){
1590 n_of_pads[1]=fTPCParam->GetNPads(isec,irow-1);
1591 n_of_pads[2]=fTPCParam->GetNPads(isec,irow);
1595 for(i1=0;i1<3;i1++){
1596 n_of_pads[i1]=fTPCParam->GetNPads(isec,irow-1+i1);
1601 // Integrated signal for this row
1602 // and a single track signal
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
1610 TMatrix &Total = *m1;
1612 // Array of pointers to the label-signal list
1614 Int_t NofDigits = n_of_pads[iFlag]*MAXTBKT; // number of digits for this row
1616 Float_t **pList = new Float_t* [NofDigits];
1620 for(lp=0;lp<NofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1623 // Straight signal and cross-talk, cross-talk is integrated over all
1624 // tracks and added to the signal at the very end
1628 for (i1=0;i1<nTracks[iFlag];i1++){
1630 m2->Zero(); // clear single track signal matrix
1632 Float_t TrackLabel =
1633 GetSignal(rowTriplet[iFlag],i1,n_of_pads[iFlag],m2,m1,IndexRange);
1635 GetList(TrackLabel,n_of_pads[iFlag],m2,IndexRange,pList);
1640 // Cross talk from the neighbouring pad-rows
1643 TMatrix *m3 = new TMatrix(0,n_of_pads[iFlag]-1,0,MAXTBKT-1); // cross-talk
1645 TMatrix &Cross = *m3;
1649 // cross-talk from the outer row only (first pad row)
1651 GetCrossTalk(0,rowTriplet[1],nTracks[1],n_of_pads,m3);
1654 else if(iFlag == 2){
1656 // cross-talk from the inner row only (last pad row)
1658 GetCrossTalk(2,rowTriplet[1],nTracks[1],n_of_pads,m3);
1663 // cross-talk from both inner and outer rows
1665 GetCrossTalk(3,rowTriplet[0],nTracks[0],n_of_pads,m3); // inner
1666 GetCrossTalk(4,rowTriplet[2],nTracks[2],n_of_pads,m3); //outer
1669 Total += Cross; // add the cross-talk
1672 // Convert analog signal to ADC counts
1679 for(Int_t ip=0;ip<n_of_pads[iFlag];ip++){
1680 for(Int_t it=0;it<MAXTBKT;it++){
1682 Float_t q = Total(ip,it);
1684 Int_t gi =it*n_of_pads[iFlag]+ip; // global index
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
1691 if(q <zerosup) continue; // do not fill zeros
1692 if(q > adc_sat) q = adc_sat; // saturation
1695 // "real" signal or electronic noise (list = -1)?
1698 for(Int_t j1=0;j1<3;j1++){
1699 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1706 digits[4]= (Int_t)q;
1710 AddDigit(tracks,digits);
1712 } // end of loop over time buckets
1713 } // end of lop over pads
1716 // This row has been digitized, delete nonused stuff
1719 for(lp=0;lp<NofDigits;lp++){
1720 if(pList[lp]) delete [] pList[lp];
1729 } // end of DigitizeRow
1730 //_____________________________________________________________________________
1731 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, Int_t np, TMatrix *m1, TMatrix *m2,
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 //---------------------------------------------------------------
1741 //-----------------------------------------------------------------
1742 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1743 //-----------------------------------------------------------------
1746 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1747 AliTPCPRF2D * fPRF2D = &(fDigParam->GetPRF2D());
1748 AliTPCRF1D * fRF = &(fDigParam->GetRF());
1750 int((z_end+6*fTPCParam->GetZSigma())/fTPCParam->GetZWidth())+1;
1752 //to make the code faster we put parameters to the stack
1754 Float_t zwidth = fTPCParam->GetZWidth();
1755 Float_t zwidthm1 =1./zwidth;
1757 tv = (TVector*)p1->At(ntr); // pointer to a track
1760 Float_t label = v(0);
1762 Int_t CentralPad = (np-1)/2;
1764 Int_t nElectrons = (tv->GetNrows()-1)/4;
1765 Float_t range=((np-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth(); // pad range
1767 range -= 0.5; // dead zone, 5mm from the edge, according to H.G. Fischer
1769 Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1772 Float_t PadSignal[7]; // signal from a single electron
1774 TMatrix &signal = *m1;
1775 TMatrix &total = *m2;
1778 IndexRange[0]=9999; // min pad
1779 IndexRange[1]=-1; // max pad
1780 IndexRange[2]=9999; //min time
1781 IndexRange[3]=-1; // max time
1784 // Loop over all electrons
1787 for(Int_t nel=0; nel<nElectrons; nel++){
1789 Float_t xwire = v(idx+1);
1790 Float_t y = v(idx+2);
1791 Float_t z = v(idx+3);
1794 Float_t absy=TMath::Abs(y);
1796 if(absy < 0.5*fTPCParam->GetPadPitchWidth()){
1797 PadNumber=CentralPad;
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);
1803 else continue; // electron out of pad-range , lost at the sector edge
1805 Float_t aval = (absy<range-0.5) ? v(idx+4):v(idx+4)*IneffFactor;
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();
1814 Int_t LeftPad = TMath::Max(0,PadNumber-3);
1815 Int_t RightPad = TMath::Min(np-1,PadNumber+3);
1817 Int_t pmin=LeftPad-PadNumber+3; // lower index of the pad_signal vector
1818 Int_t pmax=RightPad-PadNumber+3; // upper index
1820 Float_t z_drift = z*zwidthm1;
1821 Float_t z_offset = z_drift-(Int_t)z_drift;
1823 Int_t FirstBucket = (Int_t)z_drift; // numbering starts from "0"
1826 // loop over time bins (4 bins is enough - 3 sigma truncated Gaussian)
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
1834 IndexRange[2]=TMath::Min(IndexRange[2],TrueTime); // min time
1835 IndexRange[3]=TMath::Max(IndexRange[3],TrueTime); // max time
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
1848 return label; // returns track label when finished
1851 //_____________________________________________________________________________
1852 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *IndexRange,
1855 //----------------------------------------------------------------------
1856 // Updates the list of tracks contributing to digits for a given row
1857 //----------------------------------------------------------------------
1859 //-----------------------------------------------------------------
1860 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1861 //-----------------------------------------------------------------
1863 TMatrix &signal = *m;
1865 // lop over nonzero digits
1867 for(Int_t it=IndexRange[2];it<IndexRange[3]+1;it++){
1868 for(Int_t ip=IndexRange[0];ip<IndexRange[1]+1;ip++){
1871 Int_t GlobalIndex = it*np+ip; // GlobalIndex starts from 0!
1873 if(!pList[GlobalIndex]){
1876 // Create new list (6 elements - 3 signals and 3 labels),
1877 // but only if the signal is > 0.
1880 if(signal(ip,it)>0.){
1882 pList[GlobalIndex] = new Float_t [6];
1886 *pList[GlobalIndex] = -1.;
1887 *(pList[GlobalIndex]+1) = -1.;
1888 *(pList[GlobalIndex]+2) = -1.;
1889 *(pList[GlobalIndex]+3) = -1.;
1890 *(pList[GlobalIndex]+4) = -1.;
1891 *(pList[GlobalIndex]+5) = -1.;
1894 *pList[GlobalIndex] = label;
1895 *(pList[GlobalIndex]+3) = signal(ip,it);}
1899 // check the signal magnitude
1901 Float_t highest = *(pList[GlobalIndex]+3);
1902 Float_t middle = *(pList[GlobalIndex]+4);
1903 Float_t lowest = *(pList[GlobalIndex]+5);
1906 // compare the new signal with already existing list
1909 if(signal(ip,it)<lowest) continue; // neglect this track
1913 if (signal(ip,it)>highest){
1914 *(pList[GlobalIndex]+5) = middle;
1915 *(pList[GlobalIndex]+4) = highest;
1916 *(pList[GlobalIndex]+3) = signal(ip,it);
1918 *(pList[GlobalIndex]+2) = *(pList[GlobalIndex]+1);
1919 *(pList[GlobalIndex]+1) = *pList[GlobalIndex];
1920 *pList[GlobalIndex] = label;
1922 else if (signal(ip,it)>middle){
1923 *(pList[GlobalIndex]+5) = middle;
1924 *(pList[GlobalIndex]+4) = signal(ip,it);
1926 *(pList[GlobalIndex]+2) = *(pList[GlobalIndex]+1);
1927 *(pList[GlobalIndex]+1) = label;
1930 *(pList[GlobalIndex]+5) = signal(ip,it);
1931 *(pList[GlobalIndex]+2) = label;
1935 } // end of loop over pads
1936 } // end of loop over time bins
1942 //___________________________________________________________________
1943 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1944 Stat_t ntracks,TObjArray **row)
1947 //-----------------------------------------------------------------
1948 // Prepares the sector digitization, creates the vectors of
1949 // tracks for each row of this sector. The track vector
1950 // contains the track label and the position of electrons.
1951 //-----------------------------------------------------------------
1953 //-----------------------------------------------------------------
1954 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1955 //-----------------------------------------------------------------
1957 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1961 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1963 //----------------------------------------------
1964 // Create TObjArray-s, one for each row,
1965 // each TObjArray will store the TVectors
1966 // of electrons, one TVector per each track.
1967 //----------------------------------------------
1969 for(i=0; i<nrows; i++){
1970 row[i] = new TObjArray;
1972 Int_t *n_of_electrons = new Int_t [nrows]; // electron counter for each row
1973 TVector **tr = new TVector* [nrows]; //pointers to the track vectors
1975 //--------------------------------------------------------------------
1976 // Loop over tracks, the "track" contains the full history
1977 //--------------------------------------------------------------------
1979 Int_t previousTrack,currentTrack;
1980 previousTrack = -1; // nothing to store so far!
1982 for(Int_t track=0;track<ntracks;track++){
1986 TH->GetEvent(track); // get next track
1987 Int_t nhits = fHits->GetEntriesFast(); // get number of hits for this track
1989 if(nhits == 0) continue; // no hits in the TPC for this track
1991 //--------------------------------------------------------------
1993 //--------------------------------------------------------------
1995 for(Int_t hit=0;hit<nhits;hit++){
1997 tpcHit = (AliTPChit*)fHits->UncheckedAt(hit); // get a pointer to a hit
1999 Int_t sector=tpcHit->fSector; // sector number
2000 if(sector != isec) continue; //terminate iteration
2002 currentTrack = tpcHit->fTrack; // track number
2003 if(currentTrack != previousTrack){
2005 // store already filled fTrack
2007 for(i=0;i<nrows;i++){
2008 if(previousTrack != -1){
2009 if(n_of_electrons[i]>0){
2010 TVector &v = *tr[i];
2011 v(0) = previousTrack;
2012 tr[i]->ResizeTo(3*n_of_electrons[i]+1); // shrink if necessary
2016 delete tr[i]; // delete empty TVector
2021 n_of_electrons[i]=0;
2022 tr[i] = new TVector(361); // TVectors for the next fTrack
2024 } // end of loop over rows
2026 previousTrack=currentTrack; // update track label
2029 Int_t QI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2031 //---------------------------------------------------
2032 // Calculate the electron attachment probability
2033 //---------------------------------------------------
2035 Float_t time = 1.e6*(z_end-TMath::Abs(tpcHit->fZ))/fTPCParam->GetDriftV();
2037 Float_t AttProb = fTPCParam->GetAttCoef()*
2038 fTPCParam->GetOxyCont()*time; // fraction!
2040 //-----------------------------------------------
2041 // Loop over electrons
2042 //-----------------------------------------------
2044 for(Int_t nel=0;nel<QI;nel++){
2045 // skip if electron lost due to the attachment
2046 if((gRandom->Rndm(0)) < AttProb) continue; // electron lost!
2050 ElDiff(xyz); // Appply the diffusion
2052 fTPCParam->XYZtoCRXYZ(xyz,isec,row_number,3);
2054 //transform position to local coordinates
2055 //option 3 means that we calculate x position relative to
2058 if ((row_number<0)||row_number>=fTPCParam->GetNRow(isec)) continue;
2059 ExB(xyz); // ExB effect at the sense wires
2060 n_of_electrons[row_number]++;
2061 //----------------------------------
2062 // Expand vector if necessary
2063 //----------------------------------
2064 if(n_of_electrons[row_number]>120){
2065 Int_t range = tr[row_number]->GetNrows();
2066 if(n_of_electrons[row_number] > (range-1)/3){
2067 tr[row_number]->ResizeTo(range+150); // Add 50 electrons
2071 TVector &v = *tr[row_number];
2072 Int_t idx = 3*n_of_electrons[row_number]-2;
2074 v(idx)= xyz[0]; // X
2075 v(idx+1)=xyz[1]; // Y (along the pad-row)
2076 v(idx+2)=xyz[2]; // Z
2078 } // end of loop over electrons
2080 } // end of loop over hits
2081 } // end of loop over tracks
2084 // store remaining track (the last one) if not empty
2087 for(i=0;i<nrows;i++){
2088 if(n_of_electrons[i]>0){
2089 TVector &v = *tr[i];
2090 v(0) = previousTrack;
2091 tr[i]->ResizeTo(3*n_of_electrons[i]+1); // shrink if necessary
2101 delete [] n_of_electrons;
2103 } // end of MakeSector
2104 //_____________________________________________________________________________
2105 void AliTPC::GetCrossTalk (Int_t iFlag,TObjArray *p,Int_t ntracks,Int_t *npads,
2109 //-------------------------------------------------------------
2110 // Calculates the cross-talk from one row (inner or outer)
2111 //-------------------------------------------------------------
2113 //-----------------------------------------------------------------
2114 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2115 //-----------------------------------------------------------------
2118 // iFlag=2 & 3 --> cross-talk from the inner row
2119 // iFlag=0 & 4 --> cross-talk from the outer row
2121 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
2122 AliTPCPRF2D * fPRF2D = &(fDigParam->GetPRF2D());
2123 AliTPCRF1D * fRF = &(fDigParam->GetRF());
2125 int((z_end+6*fTPCParam->GetZSigma())/fTPCParam->GetZWidth())+1;
2127 //to make code faster
2129 Float_t zwidth = fTPCParam->GetZWidth();
2130 Float_t zwidthm1 =1/fTPCParam->GetZWidth();
2135 Int_t nPadsSignal; // for this pads the signal is calculated
2136 Float_t range; // sense wire range
2139 Float_t IneffFactor=0.5; // gain inefficiency close to the edges
2143 nPadsSignal = npads[1];
2144 range = ((npads[1]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2145 nPadsDiff = (npads[1]-npads[0])/2;
2147 else if (iFlag == 2){
2149 nPadsSignal = npads[2];
2150 range = ((npads[1]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2153 else if (iFlag == 3){
2155 nPadsSignal = npads[1];
2156 range = ((npads[0]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2161 nPadsSignal = npads[2];
2162 range = ((npads[2]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2163 nPadsDiff = (npads[2]-npads[1])/2;
2166 range-=0.5; // dead zone close to the edges
2169 TMatrix &signal = *m;
2171 Int_t CentralPad = (nPadsSignal-1)/2;
2172 Float_t PadSignal[7]; // signal from a single electron
2174 for(Int_t nt=0;nt<ntracks;nt++){
2175 tv=(TVector*)p->At(nt); // pointer to a track
2177 Int_t nElectrons = (tv->GetNrows()-1)/4;
2178 // Loop over electrons
2179 for(Int_t nel=0; nel<nElectrons; nel++){
2183 if (iFlag==0) xwire+=fTPCParam->GetPadPitchLength();
2184 if (iFlag==2) xwire-=fTPCParam->GetPadPitchLength();
2185 if (iFlag==3) xwire-=fTPCParam->GetPadPitchLength();
2186 if (iFlag==4) xwire+=fTPCParam->GetPadPitchLength();
2188 // electron acceptance for the cross-talk (only the closest wire)
2190 Float_t dxMax = fTPCParam->GetPadPitchLength()*0.5+fTPCParam->GetWWPitch();
2191 if(TMath::Abs(xwire)>dxMax) continue;
2193 Float_t y = v(idx+2);
2194 Float_t z = v(idx+3);
2195 Float_t absy=TMath::Abs(y);
2197 if(absy < 0.5*fTPCParam->GetPadPitchWidth()){
2198 PadNumber=CentralPad;
2200 else if (absy < range){
2201 PadNumber=(Int_t) ((absy-0.5*fTPCParam->GetPadPitchWidth())/fTPCParam->GetPadPitchWidth() +1.);
2202 PadNumber=(Int_t) (TMath::Sign((Float_t)PadNumber, y)+CentralPad);
2204 else continue; // electron out of sense wire range, lost at the sector edge
2206 Float_t aval = (absy<range-0.5) ? v(idx+4):v(idx+4)*IneffFactor;
2208 Float_t dist = y - (Float_t)(PadNumber-CentralPad)*fTPCParam->GetPadPitchWidth();
2210 for (Int_t i=0;i<7;i++){
2211 PadSignal[i]=fPRF2D->GetPRF(-dist+(i-3)*fTPCParam->GetPadPitchWidth(),xwire)*aval;
2213 PadSignal[i] *= fTPCParam->GetPadCoupling();
2217 Int_t LeftPad = TMath::Max(0,PadNumber-3);
2218 Int_t RightPad = TMath::Min(nPadsSignal-1,PadNumber+3);
2220 Int_t pmin=LeftPad-PadNumber+3; // lower index of the pad_signal vector
2221 Int_t pmax=RightPad-PadNumber+3; // upper index
2224 Float_t z_drift = z*zwidthm1;
2225 Float_t z_offset = z_drift-(Int_t)z_drift;
2227 Int_t FirstBucket = (Int_t)z_drift; // numbering starts from "0"
2229 for (Int_t i2=0;i2<4;i2++){
2230 Int_t TrueTime = FirstBucket+i2; // current time bucket
2231 Float_t dz = (Float_t(i2)+1.- z_offset)*zwidth;
2232 Float_t ampl = fRF->GetRF(dz);
2233 if((TrueTime>MAXTBKT-1)) break; // beyond the time range
2236 // loop over pads, from pmin to pmax
2238 for(Int_t i3=pmin;i3<pmax+1;i3++){
2239 Int_t TruePad = LeftPad+i3-pmin;
2241 if(TruePad<nPadsDiff || TruePad > nPadsSignal-nPadsDiff-1) continue;
2243 TruePad -= nPadsDiff;
2244 signal(TruePad,TrueTime)+=(PadSignal[i3]*ampl); // not converted to charge!
2246 } // end of loop over pads
2247 } // end of loop over time bins
2249 } // end of loop over electrons
2251 } // end of loop over tracks
2253 } // end of GetCrossTalk
2257 //_____________________________________________________________________________
2261 // Initialise TPC detector after definition of geometry
2266 for(i=0;i<35;i++) printf("*");
2267 printf(" TPC_INIT ");
2268 for(i=0;i<35;i++) printf("*");
2271 for(i=0;i<80;i++) printf("*");
2275 //_____________________________________________________________________________
2276 void AliTPC::MakeBranch(Option_t* option)
2279 // Create Tree branches for the TPC.
2281 Int_t buffersize = 4000;
2282 char branchname[10];
2283 sprintf(branchname,"%s",GetName());
2285 AliDetector::MakeBranch(option);
2287 char *D = strstr(option,"D");
2289 if (fDigits && gAlice->TreeD() && D) {
2290 gAlice->TreeD()->Branch(branchname,&fDigits, buffersize);
2291 printf("Making Branch %s for digits\n",branchname);
2294 char *R = strstr(option,"R");
2296 if (fClusters && gAlice->TreeR() && R) {
2297 gAlice->TreeR()->Branch(branchname,&fClusters, buffersize);
2298 printf("Making Branch %s for Clusters\n",branchname);
2302 //_____________________________________________________________________________
2303 void AliTPC::ResetDigits()
2306 // Reset number of digits and the digits array for this detector
2310 // if (fDigits) fDigits->Clear();
2311 if (fDigParam->GetArray()!=0) fDigParam->GetArray()->Clear();
2313 if (fClusters) fClusters->Clear();
2316 //_____________________________________________________________________________
2317 void AliTPC::SetSecAL(Int_t sec)
2319 //---------------------------------------------------
2320 // Activate/deactivate selection for lower sectors
2321 //---------------------------------------------------
2323 //-----------------------------------------------------------------
2324 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2325 //-----------------------------------------------------------------
2330 //_____________________________________________________________________________
2331 void AliTPC::SetSecAU(Int_t sec)
2333 //----------------------------------------------------
2334 // Activate/deactivate selection for upper sectors
2335 //---------------------------------------------------
2337 //-----------------------------------------------------------------
2338 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2339 //-----------------------------------------------------------------
2344 //_____________________________________________________________________________
2345 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2347 //----------------------------------------
2348 // Select active lower sectors
2349 //----------------------------------------
2351 //-----------------------------------------------------------------
2352 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2353 //-----------------------------------------------------------------
2363 //_____________________________________________________________________________
2364 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2365 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2366 Int_t s11 , Int_t s12)
2368 //--------------------------------
2369 // Select active upper sectors
2370 //--------------------------------
2372 //-----------------------------------------------------------------
2373 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2374 //-----------------------------------------------------------------
2390 //_____________________________________________________________________________
2391 void AliTPC::SetSens(Int_t sens)
2394 //-------------------------------------------------------------
2395 // Activates/deactivates the sensitive strips at the center of
2396 // the pad row -- this is for the space-point resolution calculations
2397 //-------------------------------------------------------------
2399 //-----------------------------------------------------------------
2400 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2401 //-----------------------------------------------------------------
2406 void AliTPC::SetSide(Float_t side)
2411 //____________________________________________________________________________
2412 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2413 Float_t p2,Float_t p3)
2428 //_____________________________________________________________________________
2429 void AliTPC::Streamer(TBuffer &R__b)
2432 // Stream an object of class AliTPC.
2434 if (R__b.IsReading()) {
2435 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
2436 AliDetector::Streamer(R__b);
2437 if (R__v < 2) return;
2441 fClustersIndex = new Int_t[fNsectors+1];
2442 fDigitsIndex = new Int_t[fNsectors+1];
2444 R__b.WriteVersion(AliTPC::IsA());
2445 AliDetector::Streamer(R__b);
2452 ClassImp(AliTPCcluster)
2454 //_____________________________________________________________________________
2455 AliTPCcluster::AliTPCcluster(Float_t *hits, Int_t *lab)
2458 // Creates a simulated cluster for the TPC
2460 fTracks[0] = lab[0];
2461 fTracks[1] = lab[1];
2462 fTracks[2] = lab[2];
2472 //_____________________________________________________________________________
2473 void AliTPCcluster::GetXYZ(Float_t *x, const AliTPCParam *par) const
2476 // Transformation from local to global coordinate system
2478 x[0]=par->GetPadRowRadii(fSector,fPadRow);
2481 par->CRXYZtoXYZ(x,fSector,fPadRow,1);
2485 //_____________________________________________________________________________
2486 Int_t AliTPCcluster::Compare(TObject * o)
2489 // compare two clusters according y coordinata
2491 AliTPCcluster *cl= (AliTPCcluster *)o;
2492 if (fY<cl->fY) return -1;
2493 if (fY==cl->fY) return 0;
2497 Bool_t AliTPCcluster::IsSortable() const
2500 //make AliTPCcluster sortabale
2507 ClassImp(AliTPCdigit)
2509 //_____________________________________________________________________________
2510 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2514 // Creates a TPC digit object
2516 fSector = digits[0];
2517 fPadRow = digits[1];
2520 fSignal = digits[4];
2526 //_____________________________________________________________________________
2527 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2531 // Creates a TPC hit object
2542 ClassImp(AliTPCtrack)
2544 //_____________________________________________________________________________
2545 AliTPCtrack::AliTPCtrack(Float_t *hits)
2548 // Default creator for a TPC reconstructed track object
2550 fX=hits[0]; // This is dummy code !
2553 AliTPCtrack::AliTPCtrack(const AliTPCcluster *c,const TVector& xx,
2554 const TMatrix& CC, Double_t xref, Double_t alpha):
2555 x(xx),C(CC),fClusters(200)
2557 //-----------------------------------------------------------------
2558 // This is the main track constructor.
2560 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2561 //-----------------------------------------------------------------
2565 fClusters.AddLast((AliTPCcluster*)(c));
2568 //_____________________________________________________________________________
2569 AliTPCtrack::AliTPCtrack(const AliTPCtrack& t) : x(t.x), C(t.C),
2570 fClusters(t.fClusters.GetEntriesFast())
2572 //-----------------------------------------------------------------
2573 // This is a track copy constructor.
2575 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2576 //-----------------------------------------------------------------
2580 int n=t.fClusters.GetEntriesFast();
2581 for (int i=0; i<n; i++) fClusters.AddLast(t.fClusters.UncheckedAt(i));
2584 //_____________________________________________________________________________
2585 Int_t AliTPCtrack::Compare(TObject *o) {
2586 //-----------------------------------------------------------------
2587 // This function compares tracks according to their curvature.
2589 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2590 //-----------------------------------------------------------------
2591 AliTPCtrack *t=(AliTPCtrack*)o;
2592 Double_t co=t->GetSigmaY2();
2593 Double_t c =GetSigmaY2();
2595 else if (c<co) return -1;
2599 //_____________________________________________________________________________
2600 int AliTPCtrack::PropagateTo(Double_t xk,Double_t x0,Double_t rho,Double_t pm)
2602 //-----------------------------------------------------------------
2603 // This function propagates a track to a reference plane x=xk.
2605 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2606 //-----------------------------------------------------------------
2607 if (TMath::Abs(x(2)*xk - x(3)) >= 0.999) {
2608 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Propagation failed !\n";
2612 Double_t x1=fX, x2=x1+0.5*(xk-x1), dx=x2-x1, y1=x(0), z1=x(1);
2613 Double_t c1=x(2)*x1 - x(3), r1=sqrt(1.- c1*c1);
2614 Double_t c2=x(2)*x2 - x(3), r2=sqrt(1.- c2*c2);
2616 x(0) += dx*(c1+c2)/(r1+r2);
2617 x(1) += dx*(c1+c2)/(c1*r2 + c2*r1)*x(4);
2619 TMatrix F(5,5); F.UnitMatrix();
2620 Double_t rr=r1+r2, cc=c1+c2, xx=x1+x2;
2621 F(0,2)= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
2622 F(0,3)=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
2623 Double_t cr=c1*r2+c2*r1;
2624 F(1,2)= dx*x(4)*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
2625 F(1,3)=-dx*x(4)*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
2627 TMatrix tmp(F,TMatrix::kMult,C);
2628 C.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2632 //Multiple scattering******************
2633 Double_t ey=x(2)*fX - x(3);
2634 Double_t ex=sqrt(1-ey*ey);
2636 TMatrix Q(5,5); Q=0.;
2637 Q(2,2)=ez*ez+ey*ey; Q(2,3)=-ex*ey; Q(2,4)=-ex*ez;
2638 Q(3,2)=Q(2,3); Q(3,3)= ez*ez+ex*ex; Q(3,4)=-ey*ez;
2639 Q(4,2)=Q(2,4); Q(4,3)= Q(3,4); Q(4,4)=1.;
2642 F(2,2)=-x(2)*ex; F(2,3)=-x(2)*ey;
2643 F(3,2)=-ex*(x(2)*fX-ey); F(3,3)=-(1.+ x(2)*fX*ey - ey*ey);
2644 F(4,2)=-ez*ex; F(4,3)=-ez*ey; F(4,4)=1.;
2647 Q.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2649 Double_t p2=GetPt()*GetPt()*(1.+x(4)*x(4));
2650 Double_t beta2=p2/(p2 + pm*pm);
2651 Double_t d=sqrt((x1-fX)*(x1-fX)+(y1-x(0))*(y1-x(0))+(z1-x(1))*(z1-x(1)));
2653 Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho;
2657 //Energy losses************************
2658 Double_t dE=0.153e-3/beta2*(log(5940*beta2/(1-beta2)) - beta2)*d*rho;
2659 if (x1 < x2) dE=-dE;
2660 x(2)*=(1.- sqrt(p2+pm*pm)/p2*dE);
2661 //x(3)*=(1.- sqrt(p2+pm*pm)/p2*dE);
2663 x1=fX; x2=xk; y1=x(0); z1=x(1);
2664 c1=x(2)*x1 - x(3); r1=sqrt(1.- c1*c1);
2665 c2=x(2)*x2 - x(3); r2=sqrt(1.- c2*c2);
2667 x(0) += dx*(c1+c2)/(r1+r2);
2668 x(1) += dx*(c1+c2)/(c1*r2 + c2*r1)*x(4);
2671 rr=r1+r2; cc=c1+c2; xx=x1+x2;
2672 F(0,2)= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
2673 F(0,3)=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
2675 F(1,2)= dx*x(4)*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
2676 F(1,3)=-dx*x(4)*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
2679 C.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2686 //_____________________________________________________________________________
2687 void AliTPCtrack::PropagateToVertex(Double_t x0,Double_t rho,Double_t pm)
2689 //-----------------------------------------------------------------
2690 // This function propagates tracks to the "vertex".
2692 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2693 //-----------------------------------------------------------------
2694 Double_t c=x(2)*fX - x(3);
2695 Double_t tgf=-x(3)/(x(2)*x(0) + sqrt(1-c*c));
2696 Double_t snf=tgf/sqrt(1.+ tgf*tgf);
2697 Double_t xv=(x(3)+snf)/x(2);
2698 PropagateTo(xv,x0,rho,pm);
2701 //_____________________________________________________________________________
2702 void AliTPCtrack::Update(const AliTPCcluster *c, Double_t chisq)
2704 //-----------------------------------------------------------------
2705 // This function associates a clusters with this track.
2707 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2708 //-----------------------------------------------------------------
2709 TMatrix H(2,5); H.UnitMatrix();
2710 TMatrix Ht(TMatrix::kTransposed,H);
2711 TVector m(2); m(0)=c->fY; m(1)=c->fZ;
2712 TMatrix V(2,2); V(0,0)=c->fSigmaY2; V(0,1)=0.; V(1,0)=0.; V(1,1)=c->fSigmaZ2;
2714 TMatrix tmp(H,TMatrix::kMult,C);
2715 TMatrix R(tmp,TMatrix::kMult,Ht); R+=V;
2717 Double_t det=(Double_t)R(0,0)*R(1,1) - (Double_t)R(0,1)*R(1,0);
2718 R(0,1)=R(0,0); R(0,0)=R(1,1); R(1,1)=R(0,1);
2719 R(1,0)*=-1; R(0,1)=R(1,0);
2724 TMatrix K(C,TMatrix::kMult,Ht); K*=R;
2727 x*=H; x-=m; x*=-1; x*=K; x+=savex;
2728 if (TMath::Abs(x(2)*fX-x(3)) >= 0.999) {
2729 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Filtering failed !\n";
2735 C.Mult(K,tmp); C-=saveC; C*=-1;
2737 fClusters.AddLast((AliTPCcluster*)c);
2741 //_____________________________________________________________________________
2742 int AliTPCtrack::Rotate(Double_t alpha)
2744 //-----------------------------------------------------------------
2745 // This function rotates this track.
2747 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2748 //-----------------------------------------------------------------
2751 Double_t x1=fX, y1=x(0);
2752 Double_t ca=cos(alpha), sa=sin(alpha);
2753 Double_t r1=x(2)*fX - x(3);
2756 x(0)=-x1*sa + y1*ca;
2757 x(3)=x(3)*ca + (x(2)*y1 + sqrt(1.- r1*r1))*sa;
2759 Double_t r2=x(2)*fX - x(3);
2760 if (TMath::Abs(r2) >= 0.999) {
2761 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Rotation failed !\n";
2765 Double_t y0=x(0) + sqrt(1.- r2*r2)/x(2);
2766 if ((x(0)-y0)*x(2) >= 0.) {
2767 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Rotation failed !!!\n";
2771 TMatrix F(5,5); F.UnitMatrix();
2774 F(3,2)=(y1 - r1*x1/sqrt(1.- r1*r1))*sa;
2775 F(3,3)= ca + sa*r1/sqrt(1.- r1*r1);
2776 TMatrix tmp(F,TMatrix::kMult,C);
2777 // Double_t dy2=C(0,0);
2778 C.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2779 // C(0,0)+=dy2*sa*sa*r1*r1/(1.- r1*r1);
2780 // C(1,1)+=dy2*sa*sa*x(4)*x(4)/(1.- r1*r1);
2785 //_____________________________________________________________________________
2786 void AliTPCtrack::UseClusters() const
2788 //-----------------------------------------------------------------
2789 // This function marks clusters associated with this track.
2791 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2792 //-----------------------------------------------------------------
2793 int num_of_clusters=fClusters.GetEntriesFast();
2794 for (int i=0; i<num_of_clusters; i++) {
2795 //if (i<=14) continue;
2796 AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(i);
2801 //_____________________________________________________________________________
2802 Double_t AliTPCtrack::GetPredictedChi2(const AliTPCcluster *c) const
2804 //-----------------------------------------------------------------
2805 // This function calculates a predicted chi2 increment.
2807 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2808 //-----------------------------------------------------------------
2809 TMatrix H(2,5); H.UnitMatrix();
2810 TVector m(2); m(0)=c->fY; m(1)=c->fZ;
2811 TMatrix V(2,2); V(0,0)=c->fSigmaY2; V(0,1)=0.; V(1,0)=0.; V(1,1)=c->fSigmaZ2;
2812 TVector res=x; res*=H; res-=m; //res*=-1;
2813 TMatrix tmp(H,TMatrix::kMult,C);
2814 TMatrix R(tmp,TMatrix::kMult,TMatrix(TMatrix::kTransposed,H)); R+=V;
2816 Double_t det=(Double_t)R(0,0)*R(1,1) - (Double_t)R(0,1)*R(1,0);
2817 if (TMath::Abs(det) < 1.e-10) {
2818 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Singular matrix !\n";
2821 R(0,1)=R(0,0); R(0,0)=R(1,1); R(1,1)=R(0,1);
2822 R(1,0)*=-1; R(0,1)=R(1,0);
2832 //_____________________________________________________________________________
2833 struct S { int lab; int max; };
2834 int AliTPCtrack::GetLabel(int nrows) const
2836 //-----------------------------------------------------------------
2837 // This function returns the track label. If label<0, this track is fake.
2839 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2840 //-----------------------------------------------------------------
2841 int num_of_clusters=fClusters.GetEntriesFast();
2842 S *s=new S[num_of_clusters];
2844 for (i=0; i<num_of_clusters; i++) s[i].lab=s[i].max=0;
2847 for (i=0; i<num_of_clusters; i++) {
2848 AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(i);
2849 lab=TMath::Abs(c->fTracks[0]);
2851 for (j=0; j<num_of_clusters; j++)
2852 if (s[j].lab==lab || s[j].max==0) break;
2858 for (i=0; i<num_of_clusters; i++)
2859 if (s[i].max>max) {max=s[i].max; lab=s[i].lab;}
2863 for (i=0; i<num_of_clusters; i++) {
2864 AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(i);
2865 if (TMath::Abs(c->fTracks[1]) == lab ||
2866 TMath::Abs(c->fTracks[2]) == lab ) max++;
2869 if (1.-float(max)/num_of_clusters > 0.10) return -lab;
2871 int tail=int(0.08*nrows);
2872 if (num_of_clusters < tail) return lab;
2875 for (i=1; i<=tail; i++) {
2876 AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(num_of_clusters-i);
2877 if (lab == TMath::Abs(c->fTracks[0]) ||
2878 lab == TMath::Abs(c->fTracks[1]) ||
2879 lab == TMath::Abs(c->fTracks[2])) max++;
2881 if (max < int(0.5*tail)) return -lab;
2886 //_____________________________________________________________________________
2887 void AliTPCtrack::GetPxPyPz(Double_t& px, Double_t& py, Double_t& pz) const
2889 //-----------------------------------------------------------------
2890 // This function returns reconstructed track momentum in the global system.
2892 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2893 //-----------------------------------------------------------------
2894 Double_t pt=TMath::Abs(GetPt()); // GeV/c
2895 Double_t r=x(2)*fX-x(3);
2896 Double_t y0=x(0) + sqrt(1.- r*r)/x(2);
2897 px=-pt*(x(0)-y0)*x(2); //cos(phi);
2898 py=-pt*(x(3)-fX*x(2)); //sin(phi);
2900 Double_t tmp=px*TMath::Cos(fAlpha) - py*TMath::Sin(fAlpha);
2901 py=px*TMath::Sin(fAlpha) + py*TMath::Cos(fAlpha);
2905 //_____________________________________________________________________________
2906 Double_t AliTPCtrack::GetdEdX(Double_t low, Double_t up) const {
2907 //-----------------------------------------------------------------
2908 // This funtion calculates dE/dX within the "low" and "up" cuts.
2910 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2911 //-----------------------------------------------------------------
2912 int ncl=fClusters.GetEntriesFast();
2914 Double_t *q=new Double_t[ncl];
2916 for (i=0; i<ncl; i++) {
2917 AliTPCcluster *cl=(AliTPCcluster*)(fClusters.UncheckedAt(i));
2918 // if (cl->fdEdX > 3000) continue;
2919 if (cl->fdEdX <= 0) continue;
2927 for (i=0; i<n-1; i++) {
2928 if (q[i]<=q[i+1]) continue;
2929 Double_t tmp=q[i]; q[i]=q[i+1]; q[i+1]=tmp;
2934 int nl=int(low*n), nu=int(up *n);
2936 for (i=nl; i<=nu; i++) dedx += q[i];
2941 //_________________________________________________________________________
2943 // Classes for internal tracking use
2944 //_________________________________________________________________________
2945 void AliTPCRow::InsertCluster(const AliTPCcluster* c) {
2946 //-----------------------------------------------------------------------
2947 // Insert a cluster into this pad row in accordence with its y-coordinate
2949 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2950 //-----------------------------------------------------------------------
2951 if (num_of_clusters==MAX_CLUSTER_PER_ROW) {
2952 cerr<<"AliTPCRow::InsertCluster(): Too many clusters !\n"; return;
2954 if (num_of_clusters==0) {clusters[num_of_clusters++]=c; return;}
2956 memmove(clusters+i+1 ,clusters+i,(num_of_clusters-i)*sizeof(AliTPCcluster*));
2957 clusters[i]=c; num_of_clusters++;
2960 int AliTPCRow::Find(Double_t y) const {
2961 //-----------------------------------------------------------------------
2962 // Return the index of the nearest cluster
2964 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2965 //-----------------------------------------------------------------------
2966 if (y <= clusters[0]->fY) return 0;
2967 if (y > clusters[num_of_clusters-1]->fY) return num_of_clusters;
2968 int b=0, e=num_of_clusters-1, m=(b+e)/2;
2969 for (; b<e; m=(b+e)/2) {
2970 if (y > clusters[m]->fY) b=m+1;