1 ///////////////////////////////////////////////////////////////////////////////
3 // Time Projection Chamber //
4 // This class contains the basic functions for the Time Projection Chamber //
5 // detector. Functions specific to one particular geometry are //
6 // contained in the derived classes //
10 <img src="picts/AliTPCClass.gif">
15 ///////////////////////////////////////////////////////////////////////////////
21 #include <TGeometry.h>
24 #include <TObjectTable.h>
25 #include "TParticle.h"
33 #include "AliTPCParam.h"
35 #include "AliTPCPRF2D.h"
36 #include "AliTPCRF1D.h"
42 //_____________________________________________________________________________
46 // Default constructor
57 fDigParam= new AliTPCD();
58 fDigits = fDigParam->GetArray();
61 //_____________________________________________________________________________
62 AliTPC::AliTPC(const char *name, const char *title)
63 : AliDetector(name,title)
66 // Standard constructor
70 // Initialise arrays of hits and digits
71 fHits = new TClonesArray("AliTPChit", 176);
72 // fDigits = new TClonesArray("AliTPCdigit",10000);
74 fDigParam= new AliTPCD;
75 fDigits = fDigParam->GetArray();
77 AliTPCParam *fTPCParam = &(fDigParam->GetParam());
80 // Initialise counters
84 fNsectors = fTPCParam->GetNSector();
87 fDigitsIndex = new Int_t[fNsectors+1];
88 fClustersIndex = new Int_t[fNsectors+1];
92 // Initialise color attributes
93 SetMarkerColor(kYellow);
96 //_____________________________________________________________________________
108 if (fDigitsIndex) delete [] fDigitsIndex;
109 if (fClustersIndex) delete [] fClustersIndex;
112 //_____________________________________________________________________________
113 void AliTPC::AddCluster(Float_t *hits, Int_t *tracks)
116 // Add a simulated cluster to the list
118 if(!fClusters) fClusters=new TClonesArray("AliTPCcluster",10000);
119 TClonesArray &lclusters = *fClusters;
120 new(lclusters[fNclusters++]) AliTPCcluster(hits,tracks);
123 //_____________________________________________________________________________
124 void AliTPC::AddCluster(const AliTPCcluster &c)
127 // Add a simulated cluster copy to the list
129 if(!fClusters) fClusters=new TClonesArray("AliTPCcluster",10000);
130 TClonesArray &lclusters = *fClusters;
131 new(lclusters[fNclusters++]) AliTPCcluster(c);
134 //_____________________________________________________________________________
135 void AliTPC::AddDigit(Int_t *tracks, Int_t *digits)
138 // Add a TPC digit to the list
140 // TClonesArray &ldigits = *fDigits;
142 TClonesArray &ldigits = *fDigParam->GetArray();
143 new(ldigits[fNdigits++]) AliTPCdigit(tracks,digits);
146 //_____________________________________________________________________________
147 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
150 // Add a hit to the list
152 TClonesArray &lhits = *fHits;
153 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
156 //_____________________________________________________________________________
157 void AliTPC::AddTrack(Float_t *hits)
160 // Add a track to the list of tracks
162 TClonesArray <racks = *fTracks;
163 new(ltracks[fNtracks++]) AliTPCtrack(hits);
166 //_____________________________________________________________________________
167 void AliTPC::AddTrack(const AliTPCtrack& t)
170 // Add a track copy to the list of tracks
172 if(!fTracks) fTracks=new TClonesArray("AliTPCtrack",10000);
173 TClonesArray <racks = *fTracks;
174 new(ltracks[fNtracks++]) AliTPCtrack(t);
177 //_____________________________________________________________________________
178 void AliTPC::BuildGeometry()
181 // Build TPC ROOT TNode geometry for the event display
186 const int kColorTPC=19;
187 char name[5], title[25];
188 const Double_t kDegrad=TMath::Pi()/180;
189 const Double_t kRaddeg=180./TMath::Pi();
191 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
193 Float_t InnerOpenAngle = fTPCParam->GetInnerAngle();
194 Float_t OuterOpenAngle = fTPCParam->GetOuterAngle();
196 Float_t InnerAngleShift = fTPCParam->GetInnerAngleShift();
197 Float_t OuterAngleShift = fTPCParam->GetOuterAngleShift();
199 Int_t nLo = fTPCParam->GetNInnerSector()/2;
200 Int_t nHi = fTPCParam->GetNOuterSector()/2;
202 const Double_t loAng = (Double_t)TMath::Nint(InnerOpenAngle*kRaddeg);
203 const Double_t hiAng = (Double_t)TMath::Nint(OuterOpenAngle*kRaddeg);
204 const Double_t loAngSh = (Double_t)TMath::Nint(InnerAngleShift*kRaddeg);
205 const Double_t hiAngSh = (Double_t)TMath::Nint(OuterAngleShift*kRaddeg);
208 const Double_t loCorr = 1/TMath::Cos(0.5*loAng*kDegrad);
209 const Double_t hiCorr = 1/TMath::Cos(0.5*hiAng*kDegrad);
215 // Get ALICE top node
218 Top=gAlice->GetGeometry()->GetNode("alice");
222 rl = fTPCParam->GetInSecLowEdge();
223 ru = fTPCParam->GetInSecUpEdge();
227 sprintf(name,"LS%2.2d",i);
229 sprintf(title,"TPC low sector %3d",i);
232 tubs = new TTUBS(name,title,"void",rl*loCorr,ru*loCorr,250.,
233 loAng*(i-0.5)+loAngSh,loAng*(i+0.5)+loAngSh);
234 tubs->SetNumberOfDivisions(1);
236 Node = new TNode(name,title,name,0,0,0,"");
237 Node->SetLineColor(kColorTPC);
243 rl = fTPCParam->GetOuSecLowEdge();
244 ru = fTPCParam->GetOuSecUpEdge();
247 sprintf(name,"US%2.2d",i);
249 sprintf(title,"TPC upper sector %d",i);
251 tubs = new TTUBS(name,title,"void",rl*hiCorr,ru*hiCorr,250,
252 hiAng*(i-0.5)+hiAngSh,hiAng*(i+0.5)+hiAngSh);
253 tubs->SetNumberOfDivisions(1);
255 Node = new TNode(name,title,name,0,0,0,"");
256 Node->SetLineColor(kColorTPC);
263 //_____________________________________________________________________________
264 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
267 // Calculate distance from TPC to mouse on the display
273 //_____________________________________________________________________________
274 static Double_t SigmaY2(Double_t r, Double_t tgl, Double_t pt)
277 // Parametrised error of the cluster reconstruction (pad direction)
279 pt=TMath::Abs(pt)*1000.;
282 Double_t s=a_rphi - b_rphi*r*tgl + c_rphi*x*x + d_rphi*x;
283 if (s<0.4e-3) s=0.4e-3;
284 s*=1.3; //Iouri Belikov
288 //_____________________________________________________________________________
289 static Double_t SigmaZ2(Double_t r, Double_t tgl)
292 // Parametrised error of the cluster reconstruction (drift direction)
295 Double_t s=a_z - b_z*r*tgl + c_z*tgl*tgl;
296 if (s<0.4e-3) s=0.4e-3;
297 s*=1.3; //Iouri Belikov
301 //_____________________________________________________________________________
302 inline Double_t f1(Double_t x1,Double_t y1,
303 Double_t x2,Double_t y2,
304 Double_t x3,Double_t y3)
306 //-----------------------------------------------------------------
307 // Initial approximation of the track curvature
309 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
310 //-----------------------------------------------------------------
311 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
312 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
313 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
314 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
315 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
317 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
319 return -xr*yr/sqrt(xr*xr+yr*yr);
323 //_____________________________________________________________________________
324 inline Double_t f2(Double_t x1,Double_t y1,
325 Double_t x2,Double_t y2,
326 Double_t x3,Double_t y3)
328 //-----------------------------------------------------------------
329 // Initial approximation of the track curvature times center of curvature
331 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
332 //-----------------------------------------------------------------
333 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
334 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
335 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
336 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
337 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
339 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
341 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
344 //_____________________________________________________________________________
345 inline Double_t f3(Double_t x1,Double_t y1,
346 Double_t x2,Double_t y2,
347 Double_t z1,Double_t z2)
349 //-----------------------------------------------------------------
350 // Initial approximation of the tangent of the track dip angle
352 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
353 //-----------------------------------------------------------------
354 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
357 //_____________________________________________________________________________
358 static int FindProlongation(AliTPCtrack& t, const AliTPCSector *sec,
361 //-----------------------------------------------------------------
362 // This function tries to find a track prolongation.
364 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
365 //-----------------------------------------------------------------
366 const int ROWS_TO_SKIP=int(0.5*sec->GetNRows());
367 const Float_t MAX_CHI2=12.;
368 int try_again=ROWS_TO_SKIP;
369 Double_t alpha=sec->GetAlpha();
370 int ns=int(2*TMath::Pi()/alpha+0.5);
372 for (int nr=sec->GetRowNumber(t.GetX())-1; nr>=rf; nr--) {
373 Double_t x=sec->GetX(nr), ymax=sec->GetMaxY(nr);
374 if (!t.PropagateTo(x)) return 0;
377 Double_t max_chi2=MAX_CHI2;
378 const AliTPCRow& row=sec[s][nr];
379 Double_t sy2=SigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
380 Double_t sz2=SigmaZ2(t.GetX(),t.GetTgl());
381 Double_t road=5.*sqrt(t.GetSigmaY2() + sy2), y=t.GetY(), z=t.GetZ();
384 if (t>4) cerr<<t<<" FindProlongation warning: Too broad road !\n";
389 for (int i=row.Find(y-road); i<row; i++) {
390 AliTPCcluster* c=(AliTPCcluster*)(row[i]);
391 if (c->fY > y+road) break;
392 if (c->IsUsed()) continue;
393 if ((c->fZ - z)*(c->fZ - z) > 25.*(t.GetSigmaZ2() + sz2)) continue;
394 Double_t chi2=t.GetPredictedChi2(c);
395 if (chi2 > max_chi2) continue;
401 t.Update(cl,max_chi2);
402 Double_t ll=TMath::Sqrt((1+t.GetTgl()*t.GetTgl())/
403 (1-(t.GetC()*x-t.GetEta())*(t.GetC()*x-t.GetEta())));
404 cl->fdEdX = cl->fQ/ll;
405 try_again=ROWS_TO_SKIP;
407 if (try_again==0) break;
410 if (!t.Rotate(alpha)) return 0;
411 } else if (y <-ymax) {
413 if (!t.Rotate(-alpha)) return 0;
423 //_____________________________________________________________________________
424 static void MakeSeeds(TObjArray& seeds,const AliTPCSector *sec, int max_sec,
427 //-----------------------------------------------------------------
428 // This function creates track seeds.
430 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
431 //-----------------------------------------------------------------
432 TMatrix C(5,5); TVector x(5);
433 double alpha=sec->GetAlpha(), shift=sec->GetAlphaShift();
434 double cs=cos(alpha), sn=sin(alpha);
435 for (int ns=0; ns<max_sec; ns++) {
436 int nl=sec[(ns-1+max_sec)%max_sec][i2];
438 int nu=sec[(ns+1)%max_sec][i2];
439 const AliTPCRow& r1=sec[ns][i1];
440 for (int is=0; is < r1; is++) {
441 double x1=sec->GetX(i1), y1=r1[is]->fY, z1=r1[is]->fZ;
442 for (int js=0; js < nl+nm+nu; js++) {
443 const AliTPCcluster *cl;
445 double x2=sec->GetX(i2), y2, z2, tmp;
448 ks=(ns-1+max_sec)%max_sec;
449 const AliTPCRow& r2=sec[(ns-1+max_sec)%max_sec][i2];
451 y2=cl->fY; z2=cl->fZ;
453 y2 =-x2*sn+y2*cs; x2=tmp;
457 const AliTPCRow& r2=sec[ns][i2];
459 y2=cl->fY; z2=cl->fZ;
462 const AliTPCRow& r2=sec[(ns+1)%max_sec][i2];
464 y2=cl->fY; z2=cl->fZ;
466 y2 =x2*sn+y2*cs; x2=tmp;
469 double d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
470 if (d==0.) {cerr<<"MakeSeeds warning: Straight seed !\n"; continue;}
474 x(2)=f1(x1,y1,x2,y2,0.,0.);
475 x(3)=f2(x1,y1,x2,y2,0.,0.);
476 x(4)=f3(x1,y1,x2,y2,z1,z2);
478 if (TMath::Abs(x(2)*x1-x(3)) >= 0.999) continue;
480 if (TMath::Abs(x(4)) > 1.2) continue;
482 Double_t a=asin(x(3));
483 Double_t zv=z1 - x(4)/x(2)*(a+asin(x(2)*x1-x(3)));
484 if (TMath::Abs(zv)>33.) continue;
486 TMatrix X(6,6); X=0.;
487 X(0,0)=r1[is]->fSigmaY2; X(1,1)=r1[is]->fSigmaZ2;
488 X(2,2)=cl->fSigmaY2; X(3,3)=cl->fSigmaZ2;
489 X(4,4)=3./12.; X(5,5)=3./12.;
490 TMatrix F(5,6); F.UnitMatrix();
491 Double_t sy=sqrt(X(0,0)), sz=sqrt(X(1,1));
492 F(2,0)=(f1(x1,y1+sy,x2,y2,0.,0.)-x(2))/sy;
493 F(2,2)=(f1(x1,y1,x2,y2+sy,0.,0.)-x(2))/sy;
494 F(2,4)=(f1(x1,y1,x2,y2,0.,0.+sy)-x(2))/sy;
495 F(3,0)=(f2(x1,y1+sy,x2,y2,0.,0.)-x(3))/sy;
496 F(3,2)=(f2(x1,y1,x2,y2+sy,0.,0.)-x(3))/sy;
497 F(3,4)=(f2(x1,y1,x2,y2,0.,0.+sy)-x(3))/sy;
498 F(4,0)=(f3(x1,y1+sy,x2,y2,z1,z2)-x(4))/sy;
499 F(4,1)=(f3(x1,y1,x2,y2,z1+sz,z2)-x(4))/sz;
500 F(4,2)=(f3(x1,y1,x2,y2+sy,z1,z2)-x(4))/sy;
501 F(4,3)=(f3(x1,y1,x2,y2,z1,z2+sz)-x(4))/sz;
505 TMatrix t(F,TMatrix::kMult,X);
506 C.Mult(t,TMatrix(TMatrix::kTransposed,F));
508 AliTPCtrack *track=new AliTPCtrack(r1[is], x, C, x1, ns*alpha+shift);
509 int rc=FindProlongation(*track,sec,ns,i2);
510 if (rc<0 || *track<(i1-i2)/2) delete track;
511 else seeds.AddLast(track);
517 //_____________________________________________________________________________
518 AliTPCParam *AliTPCSector::param;
519 void AliTPC::Clusters2Tracks()
521 //-----------------------------------------------------------------
522 // This is a track finder.
524 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
525 //-----------------------------------------------------------------
526 if (!fClusters) return;
528 AliTPCParam *p=&fDigParam->GetParam();
529 AliTPCSector::SetParam(p);
531 const int nis=p->GetNInnerSector()/2;
532 AliTPCSSector *ssec=new AliTPCSSector[nis];
533 int nrow_low=ssec->GetNRows();
535 const int nos=p->GetNOuterSector()/2;
536 AliTPCLSector *lsec=new AliTPCLSector[nos];
537 int nrow_up=lsec->GetNRows();
539 int ncl=fClusters->GetEntriesFast();
541 AliTPCcluster *c=(AliTPCcluster*)fClusters->UncheckedAt(ncl);
542 Int_t sec=c->fSector, row=c->fPadRow;
544 ssec[sec%nis][row].InsertCluster(c);
547 lsec[sec%nos][row].InsertCluster(c);
551 TObjArray seeds(20000);
553 int nrows=nrow_low+nrow_up;
554 int gap=int(0.125*nrows), shift=int(0.5*gap);
555 MakeSeeds(seeds, lsec, nos, nrow_up-1, nrow_up-1-gap);
556 MakeSeeds(seeds, lsec, nos, nrow_up-1-shift, nrow_up-1-shift-gap);
561 int nseed=seeds.GetEntriesFast();
563 for (int s=0; s<nseed; s++) {
564 AliTPCtrack& t=*((AliTPCtrack*)seeds.UncheckedAt(s));
565 double alpha=t.GetAlpha();
566 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
567 if (alpha < 0. ) alpha += 2.*TMath::Pi();
568 int ns=int(alpha/lsec->GetAlpha())%nos;
570 if (!FindProlongation(t,lsec,ns)) continue;
572 alpha=t.GetAlpha() + 0.5*ssec->GetAlpha() - ssec->GetAlphaShift();
573 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
574 if (alpha < 0. ) alpha += 2.*TMath::Pi();
575 ns=int(alpha/ssec->GetAlpha())%nis; //index of the inner sector needed
577 alpha=ns*ssec->GetAlpha() - t.GetAlpha();
578 if (!t.Rotate(alpha)) continue;
580 if (!FindProlongation(t,ssec,ns)) continue;
582 if (t < int(0.4*nrows)) continue;
593 //_____________________________________________________________________________
594 void AliTPC::CreateMaterials()
596 //-----------------------------------------------
597 // Create Materials for for TPC
598 //-----------------------------------------------
600 //-----------------------------------------------------------------
601 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
602 //-----------------------------------------------------------------
604 Int_t ISXFLD=gAlice->Field()->Integ();
605 Float_t SXMGMX=gAlice->Field()->Max();
607 Float_t amat[5]; // atomic numbers
608 Float_t zmat[5]; // z
609 Float_t wmat[5]; // proportions
613 // ********************* Gases *******************
615 //--------------------------------------------------------------
617 //--------------------------------------------------------------
622 Float_t a_ne = 20.18;
627 AliMaterial(20,"Ne",a_ne,z_ne,density,999.,999.);
631 Float_t a_ar = 39.948;
636 AliMaterial(21,"Ar",a_ar,z_ar,density,999.,999.);
644 //--------------------------------------------------------------
646 //--------------------------------------------------------------
663 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
665 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
680 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
682 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
697 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
699 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
701 //----------------------------------------------------------------
702 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
703 //----------------------------------------------------------------
711 Float_t a,z,rho,absl,X0,buf[1];
714 for(nc = 0;nc<fNoComp;nc++)
717 // retrive material constants
719 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
724 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
726 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? a_pure[nnc] : amol[nnc]);
727 density += fMixtProp[nc]*rho; // density of the mixture
731 // mixture proportions by weight!
733 for(nc = 0;nc<fNoComp;nc++)
736 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
738 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ? a_pure[nnc] : amol[nnc])/am;
742 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
743 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
744 AliMixture(33,"Drift gas 3",amat,zmat,density,fNoComp,wmat);
746 AliMedium(2, "Drift gas 1", 31, 0, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
747 AliMedium(3, "Drift gas 2", 32, 0, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
748 AliMedium(4, "Drift gas 3", 33, 1, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
752 AliMaterial(24, "Air", 14.61, 7.3, .001205, 30420., 67500.);
754 AliMedium(24, "Air", 24, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
756 //----------------------------------------------------------------------
758 //----------------------------------------------------------------------
762 AliMaterial(30, "Al", 26.98, 13., 2.7, 8.9, 37.2);
764 AliMedium(0, "Al",30, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
768 AliMaterial(31, "Si", 28.086, 14.,2.33, 9.36, 999.);
770 AliMedium(7, "Al",31, 0, ISXFLD, SXMGMX, 10., .1, .1, .1, .1);
789 AliMixture(32, "Mylar",amat,zmat,density,-3,wmat);
791 AliMedium(5, "Mylar",32, 0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
798 AliMaterial(33,"C normal",12.011,6.,2.265,18.8,999.);
800 AliMedium(6,"C normal",33,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
802 // G10 for inner and outr field cage
803 // G10 is 60% SiO2 + 40% epoxy, right now I use A and Z for SiO2
819 AliMixture(34,"G10 aux.",amat,zmat,density,-2,wmat);
822 gMC->Gfmate((*fIdmate)[34],namate,a,z,rho,X0,absl,buf,nbuf);
824 Float_t thickX0 = 0.0052; // field cage in X0 units
826 Float_t thick = 2.; // in cm
830 rhoFactor = X0*thickX0/thick;
831 density = rho*rhoFactor;
833 AliMaterial(35,"G10-fc",a,z,density,999.,999.);
835 AliMedium(8,"G10-fc",35,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
837 thickX0 = 0.0027; // inner vessel (eta <0.9)
839 rhoFactor = X0*thickX0/thick;
840 density = rho*rhoFactor;
842 AliMaterial(36,"G10-iv",a,z,density,999.,999.);
844 AliMedium(9,"G10-iv",36,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
848 gMC->Gfmate((*fIdmate)[33],namate,a,z,rho,X0,absl,buf,nbuf);
850 thickX0 = 0.0133; // outer vessel
852 rhoFactor = X0*thickX0/thick;
853 density = rho*rhoFactor;
856 AliMaterial(37,"C-ov",a,z,density,999.,999.);
858 AliMedium(10,"C-ov",37,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
860 thickX0=0.015; // inner vessel (cone, eta > 0.9)
862 rhoFactor = X0*thickX0/thick;
863 density = rho*rhoFactor;
865 AliMaterial(38,"C-ivc",a,z,density,999.,999.);
867 AliMedium(11,"C-ivc",38,0, ISXFLD, SXMGMX, 10., .1, .1, .001, .01);
871 AliMedium(12,"CO2",10,0, ISXFLD, SXMGMX, 10., 999.,.1, .001, .001);
877 //_____________________________________________________________________________
879 const AliTPCdigit *dig;
881 Bin() {dig=0; idx=-1;}
884 struct PreCluster : public AliTPCcluster {
885 const AliTPCdigit* summit; //pointer to the maximum digit of this precluster
886 int idx; //index in AliTPC::fClusters
887 int npeaks; //number of peaks in this precluster
888 int ndigits; //number of digits in this precluster
891 PreCluster::PreCluster() : AliTPCcluster() {npeaks=ndigits=0;}
893 //_____________________________________________________________________________
894 static void FindPreCluster(int i,int j,int maxj,Bin *bins,PreCluster &c)
896 //-----------------------------------------------------------------
897 // This function looks for "preclusters".
899 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
900 //-----------------------------------------------------------------
901 Bin& b=bins[i*maxj+j];
902 double q=(double)TMath::Abs(b.dig->fSignal);
904 if (b.idx >= 0 && b.idx != c.idx) {
909 if (q > TMath::Abs(c.summit->fSignal)) c.summit=b.dig;
918 b.dig = 0; b.idx = c.idx;
920 if (bins[(i-1)*maxj+j].dig) FindPreCluster(i-1,j,maxj,bins,c);
921 if (bins[i*maxj+(j-1)].dig) FindPreCluster(i,j-1,maxj,bins,c);
922 if (bins[(i+1)*maxj+j].dig) FindPreCluster(i+1,j,maxj,bins,c);
923 if (bins[i*maxj+(j+1)].dig) FindPreCluster(i,j+1,maxj,bins,c);
927 //_____________________________________________________________________________
928 void AliTPC::Digits2Clusters()
930 //-----------------------------------------------------------------
931 // This is a simple cluster finder.
933 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
934 //-----------------------------------------------------------------
935 AliTPCParam *par = &(fDigParam->GetParam());
937 int inp=par->GetNPads(0, par->GetNRowLow()-1);
938 int onp=par->GetNPads(par->GetNSector()-1,par->GetNRowUp() -1);
939 const int MAXY=(inp>onp) ? inp+2 : onp+2;
940 const int MAXTBKT=int((z_end+6*par->GetZSigma())/par->GetZWidth())+1;
941 const int MAXZ=MAXTBKT+2;
942 const int THRESHOLD=20;
944 TTree *t=(TTree*)gDirectory->Get("TreeD0_Param1");
945 t->GetBranch("Digits")->SetAddress(&fDigits);
946 Int_t sectors_by_rows=(Int_t)t->GetEntries();
950 for (Int_t n=0; n<sectors_by_rows; n++) {
951 if (!t->GetEvent(n)) continue;
952 Bin *bins=new Bin[MAXY*MAXZ];
953 AliTPCdigit *dig=(AliTPCdigit*)fDigits->UncheckedAt(0);
954 int sec=dig->fSector, row=dig->fPadRow;
955 int ndigits=fDigits->GetEntriesFast();
959 int nis=par->GetNInnerSector(), nos=par->GetNOuterSector();
961 npads = par->GetNPadsLow(row);
962 sign = (sec < nis/2) ? 1 : -1;
964 npads = par->GetNPadsUp(row);
965 sign = ((sec-nis) < nos/2) ? 1 : -1;
970 for (ndig=0; ndig<ndigits; ndig++) {
971 dig=(AliTPCdigit*)fDigits->UncheckedAt(ndig);
972 int i=dig->fPad+1, j=dig->fTime+1;
974 cerr<<"AliTPC::Digits2Clusters error: pad number is out of range ! ";
975 cerr<<i<<' '<<npads<<endl;
979 cerr<<"AliTPC::Digits2Clusters error: time bucket is out of range ! ";
980 cerr<<j<<' '<<MAXTBKT<<endl;
983 if (dig->fSignal >= THRESHOLD) bins[i*MAXZ+j].dig=dig;
984 if (i==1 || i==npads || j==1 || j==MAXTBKT) dig->fSignal*=-1;
990 for (i=1; i<MAXY-1; i++) {
991 for (j=1; j<MAXZ-1; j++) {
992 if (bins[i*MAXZ+j].dig == 0) continue;
993 PreCluster c; c.summit=bins[i*MAXZ+j].dig; c.idx=ncls;
994 FindPreCluster(i, j, MAXZ, bins, c);
998 double s2 = c.fSigmaY2/c.fQ - c.fY*c.fY;
999 c.fSigmaY2 = s2 + 1./12.;
1000 c.fSigmaY2 *= par->GetPadPitchWidth()*par->GetPadPitchWidth();
1001 if (s2 != 0.) c.fSigmaY2 *= 0.17;
1003 s2 = c.fSigmaZ2/c.fQ - c.fZ*c.fZ;
1004 c.fSigmaZ2 = s2 + 1./12.;
1005 c.fSigmaZ2 *= par->GetZWidth()*par->GetZWidth();
1006 if (s2 != 0.) c.fSigmaZ2 *= 0.41;
1008 c.fY = (c.fY - 0.5 - 0.5*npads)*par->GetPadPitchWidth();
1009 c.fZ = par->GetZWidth()*c.fZ;
1010 c.fZ -= 3.*par->GetZSigma(); // PASA delay
1011 c.fZ = sign*(z_end - c.fZ);
1015 c.fTracks[0]=c.summit->fTracks[0];
1016 c.fTracks[1]=c.summit->fTracks[1];
1017 c.fTracks[2]=c.summit->fTracks[2];
1019 if (c.summit->fSignal<0) {
1024 AddCluster(c); ncls++; ncl++;
1028 for (ndig=0; ndig<ndigits; ndig++) {
1029 dig=(AliTPCdigit*)fDigits->UncheckedAt(ndig);
1030 int i=dig->fPad+1, j=dig->fTime+1;
1032 cerr<<"AliTPC::Digits2Clusters error: pad number is out of range ! ";
1033 cerr<<i<<' '<<npads<<endl;
1037 cerr<<"AliTPC::Digits2Clusters error: time bucket is out of range ! ";
1038 cerr<<j<<' '<<MAXTBKT<<endl;
1041 if (TMath::Abs(dig->fSignal)>=par->GetZeroSup()) bins[i*MAXZ+j].dig=dig;
1044 for (i=1; i<MAXY-1; i++) {
1045 for (j=1; j<MAXZ-1; j++) {
1046 if (bins[i*MAXZ+j].dig == 0) continue;
1047 PreCluster c; c.summit=bins[i*MAXZ+j].dig; c.idx=ncls;
1048 FindPreCluster(i, j, MAXZ, bins, c);
1049 if (c.ndigits < 2) continue; //noise cluster
1050 if (c.npeaks>1) continue; //overlapped cluster
1054 double s2 = c.fSigmaY2/c.fQ - c.fY*c.fY;
1055 c.fSigmaY2 = s2 + 1./12.;
1056 c.fSigmaY2 *= par->GetPadPitchWidth()*par->GetPadPitchWidth();
1057 if (s2 != 0.) c.fSigmaY2 *= 0.04;
1059 s2 = c.fSigmaZ2/c.fQ - c.fZ*c.fZ;
1060 c.fSigmaZ2 = s2 + 1./12.;
1061 c.fSigmaZ2 *= par->GetZWidth()*par->GetZWidth();
1062 if (s2 != 0.) c.fSigmaZ2 *= 0.10;
1064 c.fY = (c.fY - 0.5 - 0.5*npads)*par->GetPadPitchWidth();
1065 c.fZ = par->GetZWidth()*c.fZ;
1066 c.fZ -= 3.*par->GetZSigma(); // PASA delay
1067 c.fZ = sign*(z_end - c.fZ);
1071 c.fTracks[0]=c.summit->fTracks[0];
1072 c.fTracks[1]=c.summit->fTracks[1];
1073 c.fTracks[2]=c.summit->fTracks[2];
1075 if (c.summit->fSignal<0) {
1080 if (c.npeaks==0) {AddCluster(c); ncls++; ncl++;}
1082 new ((*fClusters)[c.idx]) AliTPCcluster(c);
1087 cerr<<"sector, row, digits, clusters: "
1088 <<sec<<' '<<row<<' '<<ndigits<<' '<<ncl<<" \r";
1096 //_____________________________________________________________________________
1097 void AliTPC::ElDiff(Float_t *xyz)
1099 //--------------------------------------------------
1100 // calculates the diffusion of a single electron
1101 //--------------------------------------------------
1103 //-----------------------------------------------------------------
1104 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1105 //-----------------------------------------------------------------
1106 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1111 driftl=z_end-TMath::Abs(xyz[2]);
1113 if(driftl<0.01) driftl=0.01;
1115 // check the attachment
1117 driftl=TMath::Sqrt(driftl);
1119 // Float_t sig_t = driftl*diff_t;
1120 //Float_t sig_l = driftl*diff_l;
1121 Float_t sig_t = driftl*fTPCParam->GetDiffT();
1122 Float_t sig_l = driftl*fTPCParam->GetDiffL();
1123 xyz[0]=gRandom->Gaus(xyz[0],sig_t);
1124 xyz[1]=gRandom->Gaus(xyz[1],sig_t);
1125 xyz[2]=gRandom->Gaus(xyz[2],sig_l);
1127 if (TMath::Abs(xyz[2])>z_end){
1128 xyz[2]=TMath::Sign(z_end,z0);
1131 Float_t eps = 0.0001;
1132 xyz[2]=TMath::Sign(eps,z0);
1136 //_____________________________________________________________________________
1137 void AliTPC::Hits2Clusters()
1139 //--------------------------------------------------------
1140 // TPC simple cluster generator from hits
1141 // obtained from the TPC Fast Simulator
1142 // The point errors are taken from the parametrization
1143 //--------------------------------------------------------
1145 //-----------------------------------------------------------------
1146 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1147 //-----------------------------------------------------------------
1149 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1150 Float_t sigma_rphi,sigma_z,cl_rphi,cl_z;
1152 TParticle *particle; // pointer to a given particle
1153 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1154 TClonesArray *Particles; //pointer to the particle list
1158 Float_t pl,pt,tanth,rpad,ratio;
1161 //---------------------------------------------------------------
1162 // Get the access to the tracks
1163 //---------------------------------------------------------------
1165 TTree *TH = gAlice->TreeH();
1166 Stat_t ntracks = TH->GetEntries();
1167 Particles=gAlice->Particles();
1169 //------------------------------------------------------------
1170 // Loop over all sectors (72 sectors)
1171 // Sectors 0-35 are lower sectors, 0-17 z>0, 18-35 z<0
1172 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1174 // First cluster for sector 0 starts at "0"
1175 //------------------------------------------------------------
1178 //fClustersIndex[0] = 0;
1181 int Nsectors=fDigParam->GetParam().GetNSector();
1182 for(Int_t isec=0; isec<Nsectors; isec++){
1184 fTPCParam->AdjustAngles(isec,cph,sph);
1186 //------------------------------------------------------------
1188 //------------------------------------------------------------
1190 for(Int_t track=0;track<ntracks;track++){
1192 TH->GetEvent(track);
1194 // Get number of the TPC hits and a pointer
1197 nhits=fHits->GetEntriesFast();
1201 for(Int_t hit=0;hit<nhits;hit++){
1202 tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
1203 sector=tpcHit->fSector; // sector number
1204 if(sector != isec) continue; //terminate iteration
1205 ipart=tpcHit->fTrack;
1206 particle=(TParticle*)Particles->UncheckedAt(ipart);
1209 if(pt < 1.e-9) pt=1.e-9;
1211 tanth = TMath::Abs(tanth);
1212 rpad=TMath::Sqrt(tpcHit->fX*tpcHit->fX + tpcHit->fY*tpcHit->fY);
1213 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
1215 // space-point resolutions
1217 sigma_rphi=SigmaY2(rpad,tanth,pt);
1218 sigma_z =SigmaZ2(rpad,tanth );
1222 cl_rphi=ac_rphi-bc_rphi*rpad*tanth+cc_rphi*ratio*ratio;
1223 cl_z=ac_z-bc_z*rpad*tanth+cc_z*tanth*tanth;
1225 // temporary protection
1227 if(sigma_rphi < 0.) sigma_rphi=0.4e-3;
1228 if(sigma_z < 0.) sigma_z=0.4e-3;
1229 if(cl_rphi < 0.) cl_rphi=2.5e-3;
1230 if(cl_z < 0.) cl_z=2.5e-5;
1233 // smearing --> rotate sectors firstly,
1234 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1236 Float_t xprim= tpcHit->fX*cph + tpcHit->fY*sph;
1237 Float_t yprim=-tpcHit->fX*sph + tpcHit->fY*cph;
1238 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigma_rphi)); // y
1239 Double_t alpha=(sector < fTPCParam->GetNInnerSector()) ?
1240 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1241 if (TMath::Abs(xyz[0]/xprim) > TMath::Tan(0.5*alpha)) xyz[0]=yprim;
1242 xyz[1]=gRandom->Gaus(tpcHit->fZ,TMath::Sqrt(sigma_z)); // z
1243 if (TMath::Abs(xyz[1]) > 250) xyz[1]=tpcHit->fZ;
1244 xyz[2]=tpcHit->fQ+1;// q; let it be not equal to zero (Y.Belikov)
1245 xyz[3]=sigma_rphi; // fSigmaY2
1246 xyz[4]=sigma_z; // fSigmaZ2
1248 // and finally add the cluster
1249 Int_t tracks[5]={tpcHit->fTrack, -1, -1, sector, tpcHit->fPadRow};
1250 AddCluster(xyz,tracks);
1252 } // end of loop over hits
1253 } // end of loop over tracks
1255 //fClustersIndex[isec] = fNclusters; // update clusters index
1257 } // end of loop over sectors
1259 //fClustersIndex[fNsectors]--; // set end of the clusters buffer
1264 void AliTPC::Hits2Digits()
1267 //----------------------------------------------------
1268 // Loop over all sectors (72 sectors)
1269 // Sectors 0-35 are lower sectors, 0-17 z>0, 18-35 z<0
1270 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1272 int Nsectors=fDigParam->GetParam().GetNSector();
1273 for(Int_t isec=0;isec<Nsectors;isec++) Hits2DigitsSector(isec);
1277 //_____________________________________________________________________________
1278 void AliTPC::Hits2DigitsSector(Int_t isec)
1280 //-------------------------------------------------------------------
1281 // TPC conversion from hits to digits.
1282 //-------------------------------------------------------------------
1284 //-----------------------------------------------------------------
1285 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1286 //-----------------------------------------------------------------
1288 //-------------------------------------------------------
1289 // Get the access to the track hits
1290 //-------------------------------------------------------
1292 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1293 TTree *TH = gAlice->TreeH(); // pointer to the hits tree
1294 Stat_t ntracks = TH->GetEntries();
1298 //-------------------------------------------
1299 // Only if there are any tracks...
1300 //-------------------------------------------
1303 // TObjArrays for three neighouring pad-rows
1305 TObjArray **rowTriplet = new TObjArray* [3];
1307 // TObjArray-s for each pad-row
1311 for (Int_t trip=0;trip<3;trip++){
1312 rowTriplet[trip]=new TObjArray;
1317 printf("*** Processing sector number %d ***\n",isec);
1319 Int_t nrows =fTPCParam->GetNRow(isec);
1321 row= new TObjArray* [nrows];
1323 MakeSector(isec,nrows,TH,ntracks,row);
1325 //--------------------------------------------------------
1326 // Digitize this sector, row by row
1327 // row[i] is the pointer to the TObjArray of TVectors,
1328 // each one containing electrons accepted on this
1329 // row, assigned into tracks
1330 //--------------------------------------------------------
1334 for (i=0;i<nrows;i++){
1336 // Triplets for i = 0 and i=1 are identical!
1337 // The same for i = nrows-1 and nrows!
1339 if(i != 1 && i != nrows-1){
1340 MakeTriplet(i,rowTriplet,row);
1343 DigitizeRow(i,isec,rowTriplet);
1347 Int_t ndig=fDigParam->GetArray()->GetEntriesFast();
1349 printf("*** Sector, row, digits %d %d %d ***\n",isec,i,ndig);
1351 ResetDigits(); // reset digits for this row after storing them
1353 } // end of the sector digitization
1355 // delete the last triplet
1357 for (i=0;i<3;i++) rowTriplet[i]->Delete();
1359 delete [] row; // delete the array of pointers to TObjArray-s
1362 } // end of Hits2Digits
1363 //_____________________________________________________________________________
1364 void AliTPC::MakeTriplet(Int_t row,
1365 TObjArray **rowTriplet, TObjArray **prow)
1367 //------------------------------------------------------------------
1368 // Makes the "triplet" of the neighbouring pad-row for the
1369 // digitization including the cross-talk between the pad-rows
1370 //------------------------------------------------------------------
1372 //-----------------------------------------------------------------
1373 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1374 //-----------------------------------------------------------------
1376 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1377 Float_t gasgain = fTPCParam->GetGasGain();
1380 Int_t nElements,nElectrons;
1385 //-------------------------------------------------------------------
1386 // pv is an "old" track, i.e. label + triplets of (x,y,z)
1387 // for each electron
1389 //-------------------------------------------------------------------
1395 if(row == 0 || row == 1){
1397 // create entire triplet for the first AND the second row
1399 nTracks[0] = prow[0]->GetEntries();
1400 nTracks[1] = prow[1]->GetEntries();
1401 nTracks[2] = prow[2]->GetEntries();
1403 for(i1=0;i1<3;i1++){
1404 nt = nTracks[i1]; // number of tracks for this row
1406 for(i2=0;i2<nt;i2++){
1407 pv = (TVector*)prow[i1]->At(i2);
1409 nElements = pv->GetNrows();
1410 nElectrons = (nElements-1)/3;
1412 tv = new TVector(4*nElectrons+1); // create TVector for a modified track
1414 v2(0)=v1(0); //track label
1416 for(nel=0;nel<nElectrons;nel++){
1419 // Avalanche, including fluctuations
1420 Int_t aval = (Int_t) (-gasgain*TMath::Log(gRandom->Rndm()));
1421 v2(idx2+1)= v1(idx1+1);
1422 v2(idx2+2)= v1(idx1+2);
1423 v2(idx2+3)= v1(idx1+3);
1424 v2(idx2+4)= (Float_t)aval; // in number of electrons!
1425 } // end of loop over electrons
1427 // Add this track to a row
1430 rowTriplet[i1]->Add(tv);
1433 } // end of loop over tracks for this row
1435 prow[i1]->Delete(); // remove "old" tracks
1436 delete prow[i1]; // delete TObjArray for this row
1437 prow[i1]=0; // set pointer to NULL
1439 } // end of loop over row triplets
1445 rowTriplet[0]->Delete(); // remove old lower row
1447 nTracks[0]=rowTriplet[1]->GetEntries(); // previous middle row
1448 nTracks[1]=rowTriplet[2]->GetEntries(); // previous upper row
1449 nTracks[2]=prow[row+1]->GetEntries(); // next row
1452 //-------------------------------------------
1453 // shift new tracks downwards
1454 //-------------------------------------------
1456 for(i1=0;i1<nTracks[0];i1++){
1457 pv=(TVector*)rowTriplet[1]->At(i1);
1458 rowTriplet[0]->Add(pv);
1461 rowTriplet[1]->Clear(); // leave tracks on the heap!!!
1463 for(i1=0;i1<nTracks[1];i1++){
1464 pv=(TVector*)rowTriplet[2]->At(i1);
1465 rowTriplet[1]->Add(pv);
1468 rowTriplet[2]->Clear(); // leave tracks on the heap!!!
1470 //---------------------------------------------
1471 // Create new upper row
1472 //---------------------------------------------
1476 for(i1=0;i1<nTracks[2];i1++){
1477 pv = (TVector*)prow[row+1]->At(i1);
1479 nElements = pv->GetNrows();
1480 nElectrons = (nElements-1)/3;
1482 tv = new TVector(4*nElectrons+1); // create TVector for a modified track
1484 v2(0)=v1(0); //track label
1486 for(nel=0;nel<nElectrons;nel++){
1490 // Avalanche, including fluctuations
1491 Int_t aval = (Int_t)
1492 (-gasgain*TMath::Log(gRandom->Rndm()));
1494 v2(idx2+1)= v1(idx1+1);
1495 v2(idx2+2)= v1(idx1+2);
1496 v2(idx2+3)= v1(idx1+3);
1497 v2(idx2+4)= (Float_t)aval; // in number of electrons!
1498 } // end of loop over electrons
1500 rowTriplet[2]->Add(tv);
1502 } // end of loop over tracks
1504 prow[row+1]->Delete(); // delete tracks for this row
1505 delete prow[row+1]; // delete TObjArray for this row
1506 prow[row+1]=0; // set a pointer to NULL
1510 } // end of MakeTriplet
1511 //_____________________________________________________________________________
1512 void AliTPC::ExB(Float_t *xyz)
1514 //------------------------------------------------
1515 // ExB at the wires and wire number calulation
1516 //------------------------------------------------
1518 //-----------------------------------------------------------------
1519 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1520 //-----------------------------------------------------------------
1521 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1524 fTPCParam->GetWire(x1); //calculate nearest wire position
1525 Float_t dx=xyz[0]-x1;
1526 xyz[1]+=dx*fTPCParam->GetOmegaTau();
1529 //_____________________________________________________________________________
1530 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rowTriplet)
1532 //-----------------------------------------------------------
1533 // Single row digitization, coupling from the neighbouring
1534 // rows taken into account
1535 //-----------------------------------------------------------
1537 //-----------------------------------------------------------------
1538 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1539 //-----------------------------------------------------------------
1542 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1543 Float_t chipgain= fTPCParam->GetChipGain();
1544 Float_t zerosup = fTPCParam->GetZeroSup();
1545 Int_t nrows =fTPCParam->GetNRow(isec);
1547 int((z_end+6*fTPCParam->GetZSigma())/fTPCParam->GetZWidth())+1;
1551 Int_t IndexRange[4];
1556 // iFlag = 0 -> inner row, iFlag = 1 -> the middle one, iFlag = 2 -> the outer one
1559 nTracks[0]=rowTriplet[0]->GetEntries(); // lower row
1560 nTracks[1]=rowTriplet[1]->GetEntries(); // middle row
1561 nTracks[2]=rowTriplet[2]->GetEntries(); // upper row
1566 n_of_pads[0]=fTPCParam->GetNPads(isec,0);
1567 n_of_pads[1]=fTPCParam->GetNPads(isec,1);
1569 else if(irow == nrows-1){
1571 n_of_pads[1]=fTPCParam->GetNPads(isec,irow-1);
1572 n_of_pads[2]=fTPCParam->GetNPads(isec,irow);
1576 for(i1=0;i1<3;i1++){
1577 n_of_pads[i1]=fTPCParam->GetNPads(isec,irow-1+i1);
1582 // Integrated signal for this row
1583 // and a single track signal
1586 TMatrix *m1 = new TMatrix(0,n_of_pads[iFlag]-1,0,MAXTBKT-1); // integrated
1587 TMatrix *m2 = new TMatrix(0,n_of_pads[iFlag]-1,0,MAXTBKT-1); // single
1591 TMatrix &Total = *m1;
1593 // Array of pointers to the label-signal list
1595 Int_t NofDigits = n_of_pads[iFlag]*MAXTBKT; // number of digits for this row
1597 Float_t **pList = new Float_t* [NofDigits];
1601 for(lp=0;lp<NofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1604 // Straight signal and cross-talk, cross-talk is integrated over all
1605 // tracks and added to the signal at the very end
1609 for (i1=0;i1<nTracks[iFlag];i1++){
1611 m2->Zero(); // clear single track signal matrix
1613 Float_t TrackLabel =
1614 GetSignal(rowTriplet[iFlag],i1,n_of_pads[iFlag],m2,m1,IndexRange);
1616 GetList(TrackLabel,n_of_pads[iFlag],m2,IndexRange,pList);
1621 // Cross talk from the neighbouring pad-rows
1624 TMatrix *m3 = new TMatrix(0,n_of_pads[iFlag]-1,0,MAXTBKT-1); // cross-talk
1626 TMatrix &Cross = *m3;
1630 // cross-talk from the outer row only (first pad row)
1632 GetCrossTalk(0,rowTriplet[1],nTracks[1],n_of_pads,m3);
1635 else if(iFlag == 2){
1637 // cross-talk from the inner row only (last pad row)
1639 GetCrossTalk(2,rowTriplet[1],nTracks[1],n_of_pads,m3);
1644 // cross-talk from both inner and outer rows
1646 GetCrossTalk(3,rowTriplet[0],nTracks[0],n_of_pads,m3); // inner
1647 GetCrossTalk(4,rowTriplet[2],nTracks[2],n_of_pads,m3); //outer
1650 Total += Cross; // add the cross-talk
1653 // Convert analog signal to ADC counts
1660 for(Int_t ip=0;ip<n_of_pads[iFlag];ip++){
1661 for(Int_t it=0;it<MAXTBKT;it++){
1663 Float_t q = Total(ip,it);
1665 Int_t gi =it*n_of_pads[iFlag]+ip; // global index
1667 q = gRandom->Gaus(q,fTPCParam->GetNoise()); // apply noise
1668 q *= (q_el*1.e15); // convert to fC
1669 q *= chipgain; // convert to mV
1670 q *= (adc_sat/dyn_range); // convert to ADC counts
1672 if(q <zerosup) continue; // do not fill zeros
1673 if(q > adc_sat) q = adc_sat; // saturation
1676 // "real" signal or electronic noise (list = -1)?
1679 for(Int_t j1=0;j1<3;j1++){
1680 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1687 digits[4]= (Int_t)q;
1691 AddDigit(tracks,digits);
1693 } // end of loop over time buckets
1694 } // end of lop over pads
1697 // This row has been digitized, delete nonused stuff
1700 for(lp=0;lp<NofDigits;lp++){
1701 if(pList[lp]) delete [] pList[lp];
1710 } // end of DigitizeRow
1711 //_____________________________________________________________________________
1712 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, Int_t np, TMatrix *m1, TMatrix *m2,
1716 //---------------------------------------------------------------
1717 // Calculates 2-D signal (pad,time) for a single track,
1718 // returns a pointer to the signal matrix and the track label
1719 // No digitization is performed at this level!!!
1720 //---------------------------------------------------------------
1722 //-----------------------------------------------------------------
1723 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1724 //-----------------------------------------------------------------
1727 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1728 AliTPCPRF2D * fPRF2D = &(fDigParam->GetPRF2D());
1729 AliTPCRF1D * fRF = &(fDigParam->GetRF());
1731 int((z_end+6*fTPCParam->GetZSigma())/fTPCParam->GetZWidth())+1;
1733 //to make the code faster we put parameters to the stack
1735 Float_t zwidth = fTPCParam->GetZWidth();
1736 Float_t zwidthm1 =1./zwidth;
1738 tv = (TVector*)p1->At(ntr); // pointer to a track
1741 Float_t label = v(0);
1743 Int_t CentralPad = (np-1)/2;
1745 Int_t nElectrons = (tv->GetNrows()-1)/4;
1746 Float_t range=((np-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth(); // pad range
1748 range -= 0.5; // dead zone, 5mm from the edge, according to H.G. Fischer
1750 Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1753 Float_t PadSignal[7]; // signal from a single electron
1755 TMatrix &signal = *m1;
1756 TMatrix &total = *m2;
1759 IndexRange[0]=9999; // min pad
1760 IndexRange[1]=-1; // max pad
1761 IndexRange[2]=9999; //min time
1762 IndexRange[3]=-1; // max time
1765 // Loop over all electrons
1768 for(Int_t nel=0; nel<nElectrons; nel++){
1770 Float_t xwire = v(idx+1);
1771 Float_t y = v(idx+2);
1772 Float_t z = v(idx+3);
1775 Float_t absy=TMath::Abs(y);
1777 if(absy < 0.5*fTPCParam->GetPadPitchWidth()){
1778 PadNumber=CentralPad;
1780 else if (absy < range){
1781 PadNumber=(Int_t) ((absy-0.5*fTPCParam->GetPadPitchWidth())/fTPCParam->GetPadPitchWidth()+1.);
1782 PadNumber=(Int_t) (TMath::Sign((Float_t)PadNumber, y)+CentralPad);
1784 else continue; // electron out of pad-range , lost at the sector edge
1786 Float_t aval = (absy<range-0.5) ? v(idx+4):v(idx+4)*IneffFactor;
1789 Float_t dist = y - (Float_t)(PadNumber-CentralPad)*fTPCParam->GetPadPitchWidth();
1790 for (Int_t i=0;i<7;i++){
1791 PadSignal[i]=fPRF2D->GetPRF(-dist+(i-3)*fTPCParam->GetPadPitchWidth(),xwire)*aval;
1792 PadSignal[i] *= fTPCParam->GetPadCoupling();
1795 Int_t LeftPad = TMath::Max(0,PadNumber-3);
1796 Int_t RightPad = TMath::Min(np-1,PadNumber+3);
1798 Int_t pmin=LeftPad-PadNumber+3; // lower index of the pad_signal vector
1799 Int_t pmax=RightPad-PadNumber+3; // upper index
1801 Float_t z_drift = z*zwidthm1;
1802 Float_t z_offset = z_drift-(Int_t)z_drift;
1804 Int_t FirstBucket = (Int_t)z_drift; // numbering starts from "0"
1807 // loop over time bins (4 bins is enough - 3 sigma truncated Gaussian)
1809 for (Int_t i2=0;i2<4;i2++){
1810 Int_t TrueTime = FirstBucket+i2; // current time bucket
1811 Float_t dz = (Float_t(i2)+1.-z_offset)*zwidth;
1812 Float_t ampl = fRF->GetRF(dz);
1813 if( (TrueTime>MAXTBKT-1) ) break; // beyond the time range
1815 IndexRange[2]=TMath::Min(IndexRange[2],TrueTime); // min time
1816 IndexRange[3]=TMath::Max(IndexRange[3],TrueTime); // max time
1818 // loop over pads, from pmin to pmax
1819 for(Int_t i3=pmin;i3<=pmax;i3++){
1820 Int_t TruePad = LeftPad+i3-pmin;
1821 IndexRange[0]=TMath::Min(IndexRange[0],TruePad); // min pad
1822 IndexRange[1]=TMath::Max(IndexRange[1],TruePad); // max pad
1823 signal(TruePad,TrueTime)+=(PadSignal[i3]*ampl); // not converted to charge!!!
1824 total(TruePad,TrueTime)+=(PadSignal[i3]*ampl); // not converted to charge!!!
1825 } // end of pads loop
1826 } // end of time loop
1827 } // end of loop over electrons
1829 return label; // returns track label when finished
1832 //_____________________________________________________________________________
1833 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *IndexRange,
1836 //----------------------------------------------------------------------
1837 // Updates the list of tracks contributing to digits for a given row
1838 //----------------------------------------------------------------------
1840 //-----------------------------------------------------------------
1841 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1842 //-----------------------------------------------------------------
1844 TMatrix &signal = *m;
1846 // lop over nonzero digits
1848 for(Int_t it=IndexRange[2];it<IndexRange[3]+1;it++){
1849 for(Int_t ip=IndexRange[0];ip<IndexRange[1]+1;ip++){
1852 Int_t GlobalIndex = it*np+ip; // GlobalIndex starts from 0!
1854 if(!pList[GlobalIndex]){
1857 // Create new list (6 elements - 3 signals and 3 labels),
1858 // but only if the signal is > 0.
1861 if(signal(ip,it)>0.){
1863 pList[GlobalIndex] = new Float_t [6];
1867 *pList[GlobalIndex] = -1.;
1868 *(pList[GlobalIndex]+1) = -1.;
1869 *(pList[GlobalIndex]+2) = -1.;
1870 *(pList[GlobalIndex]+3) = -1.;
1871 *(pList[GlobalIndex]+4) = -1.;
1872 *(pList[GlobalIndex]+5) = -1.;
1875 *pList[GlobalIndex] = label;
1876 *(pList[GlobalIndex]+3) = signal(ip,it);}
1880 // check the signal magnitude
1882 Float_t highest = *(pList[GlobalIndex]+3);
1883 Float_t middle = *(pList[GlobalIndex]+4);
1884 Float_t lowest = *(pList[GlobalIndex]+5);
1887 // compare the new signal with already existing list
1890 if(signal(ip,it)<lowest) continue; // neglect this track
1894 if (signal(ip,it)>highest){
1895 *(pList[GlobalIndex]+5) = middle;
1896 *(pList[GlobalIndex]+4) = highest;
1897 *(pList[GlobalIndex]+3) = signal(ip,it);
1899 *(pList[GlobalIndex]+2) = *(pList[GlobalIndex]+1);
1900 *(pList[GlobalIndex]+1) = *pList[GlobalIndex];
1901 *pList[GlobalIndex] = label;
1903 else if (signal(ip,it)>middle){
1904 *(pList[GlobalIndex]+5) = middle;
1905 *(pList[GlobalIndex]+4) = signal(ip,it);
1907 *(pList[GlobalIndex]+2) = *(pList[GlobalIndex]+1);
1908 *(pList[GlobalIndex]+1) = label;
1911 *(pList[GlobalIndex]+5) = signal(ip,it);
1912 *(pList[GlobalIndex]+2) = label;
1916 } // end of loop over pads
1917 } // end of loop over time bins
1923 //___________________________________________________________________
1924 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1925 Stat_t ntracks,TObjArray **row)
1928 //-----------------------------------------------------------------
1929 // Prepares the sector digitization, creates the vectors of
1930 // tracks for each row of this sector. The track vector
1931 // contains the track label and the position of electrons.
1932 //-----------------------------------------------------------------
1934 //-----------------------------------------------------------------
1935 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1936 //-----------------------------------------------------------------
1938 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
1942 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1944 //----------------------------------------------
1945 // Create TObjArray-s, one for each row,
1946 // each TObjArray will store the TVectors
1947 // of electrons, one TVector per each track.
1948 //----------------------------------------------
1950 for(i=0; i<nrows; i++){
1951 row[i] = new TObjArray;
1953 Int_t *n_of_electrons = new Int_t [nrows]; // electron counter for each row
1954 TVector **tr = new TVector* [nrows]; //pointers to the track vectors
1956 //--------------------------------------------------------------------
1957 // Loop over tracks, the "track" contains the full history
1958 //--------------------------------------------------------------------
1960 Int_t previousTrack,currentTrack;
1961 previousTrack = -1; // nothing to store so far!
1963 for(Int_t track=0;track<ntracks;track++){
1967 TH->GetEvent(track); // get next track
1968 Int_t nhits = fHits->GetEntriesFast(); // get number of hits for this track
1970 if(nhits == 0) continue; // no hits in the TPC for this track
1972 //--------------------------------------------------------------
1974 //--------------------------------------------------------------
1976 for(Int_t hit=0;hit<nhits;hit++){
1978 tpcHit = (AliTPChit*)fHits->UncheckedAt(hit); // get a pointer to a hit
1980 Int_t sector=tpcHit->fSector; // sector number
1981 if(sector != isec) continue; //terminate iteration
1983 currentTrack = tpcHit->fTrack; // track number
1984 if(currentTrack != previousTrack){
1986 // store already filled fTrack
1988 for(i=0;i<nrows;i++){
1989 if(previousTrack != -1){
1990 if(n_of_electrons[i]>0){
1991 TVector &v = *tr[i];
1992 v(0) = previousTrack;
1993 tr[i]->ResizeTo(3*n_of_electrons[i]+1); // shrink if necessary
1997 delete tr[i]; // delete empty TVector
2002 n_of_electrons[i]=0;
2003 tr[i] = new TVector(361); // TVectors for the next fTrack
2005 } // end of loop over rows
2007 previousTrack=currentTrack; // update track label
2010 Int_t QI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2012 //---------------------------------------------------
2013 // Calculate the electron attachment probability
2014 //---------------------------------------------------
2016 Float_t time = 1.e6*(z_end-TMath::Abs(tpcHit->fZ))/fTPCParam->GetDriftV();
2018 Float_t AttProb = fTPCParam->GetAttCoef()*
2019 fTPCParam->GetOxyCont()*time; // fraction!
2021 //-----------------------------------------------
2022 // Loop over electrons
2023 //-----------------------------------------------
2025 for(Int_t nel=0;nel<QI;nel++){
2026 // skip if electron lost due to the attachment
2027 if((gRandom->Rndm(0)) < AttProb) continue; // electron lost!
2031 ElDiff(xyz); // Appply the diffusion
2033 fTPCParam->XYZtoCRXYZ(xyz,isec,row_number,3);
2035 //transform position to local coordinates
2036 //option 3 means that we calculate x position relative to
2039 if ((row_number<0)||row_number>=fTPCParam->GetNRow(isec)) continue;
2040 ExB(xyz); // ExB effect at the sense wires
2041 n_of_electrons[row_number]++;
2042 //----------------------------------
2043 // Expand vector if necessary
2044 //----------------------------------
2045 if(n_of_electrons[row_number]>120){
2046 Int_t range = tr[row_number]->GetNrows();
2047 if(n_of_electrons[row_number] > (range-1)/3){
2048 tr[row_number]->ResizeTo(range+150); // Add 50 electrons
2052 TVector &v = *tr[row_number];
2053 Int_t idx = 3*n_of_electrons[row_number]-2;
2055 v(idx)= xyz[0]; // X
2056 v(idx+1)=xyz[1]; // Y (along the pad-row)
2057 v(idx+2)=xyz[2]; // Z
2059 } // end of loop over electrons
2061 } // end of loop over hits
2062 } // end of loop over tracks
2065 // store remaining track (the last one) if not empty
2068 for(i=0;i<nrows;i++){
2069 if(n_of_electrons[i]>0){
2070 TVector &v = *tr[i];
2071 v(0) = previousTrack;
2072 tr[i]->ResizeTo(3*n_of_electrons[i]+1); // shrink if necessary
2082 delete [] n_of_electrons;
2084 } // end of MakeSector
2085 //_____________________________________________________________________________
2086 void AliTPC::GetCrossTalk (Int_t iFlag,TObjArray *p,Int_t ntracks,Int_t *npads,
2090 //-------------------------------------------------------------
2091 // Calculates the cross-talk from one row (inner or outer)
2092 //-------------------------------------------------------------
2094 //-----------------------------------------------------------------
2095 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2096 //-----------------------------------------------------------------
2099 // iFlag=2 & 3 --> cross-talk from the inner row
2100 // iFlag=0 & 4 --> cross-talk from the outer row
2102 AliTPCParam * fTPCParam = &(fDigParam->GetParam());
2103 AliTPCPRF2D * fPRF2D = &(fDigParam->GetPRF2D());
2104 AliTPCRF1D * fRF = &(fDigParam->GetRF());
2106 int((z_end+6*fTPCParam->GetZSigma())/fTPCParam->GetZWidth())+1;
2108 //to make code faster
2110 Float_t zwidth = fTPCParam->GetZWidth();
2111 Float_t zwidthm1 =1/fTPCParam->GetZWidth();
2116 Int_t nPadsSignal; // for this pads the signal is calculated
2117 Float_t range; // sense wire range
2120 Float_t IneffFactor=0.5; // gain inefficiency close to the edges
2124 nPadsSignal = npads[1];
2125 range = ((npads[1]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2126 nPadsDiff = (npads[1]-npads[0])/2;
2128 else if (iFlag == 2){
2130 nPadsSignal = npads[2];
2131 range = ((npads[1]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2134 else if (iFlag == 3){
2136 nPadsSignal = npads[1];
2137 range = ((npads[0]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2142 nPadsSignal = npads[2];
2143 range = ((npads[2]-1)/2 + 0.5)*fTPCParam->GetPadPitchWidth();
2144 nPadsDiff = (npads[2]-npads[1])/2;
2147 range-=0.5; // dead zone close to the edges
2150 TMatrix &signal = *m;
2152 Int_t CentralPad = (nPadsSignal-1)/2;
2153 Float_t PadSignal[7]; // signal from a single electron
2155 for(Int_t nt=0;nt<ntracks;nt++){
2156 tv=(TVector*)p->At(nt); // pointer to a track
2158 Int_t nElectrons = (tv->GetNrows()-1)/4;
2159 // Loop over electrons
2160 for(Int_t nel=0; nel<nElectrons; nel++){
2164 if (iFlag==0) xwire+=fTPCParam->GetPadPitchLength();
2165 if (iFlag==2) xwire-=fTPCParam->GetPadPitchLength();
2166 if (iFlag==3) xwire-=fTPCParam->GetPadPitchLength();
2167 if (iFlag==4) xwire+=fTPCParam->GetPadPitchLength();
2169 // electron acceptance for the cross-talk (only the closest wire)
2171 Float_t dxMax = fTPCParam->GetPadPitchLength()*0.5+fTPCParam->GetWWPitch();
2172 if(TMath::Abs(xwire)>dxMax) continue;
2174 Float_t y = v(idx+2);
2175 Float_t z = v(idx+3);
2176 Float_t absy=TMath::Abs(y);
2178 if(absy < 0.5*fTPCParam->GetPadPitchWidth()){
2179 PadNumber=CentralPad;
2181 else if (absy < range){
2182 PadNumber=(Int_t) ((absy-0.5*fTPCParam->GetPadPitchWidth())/fTPCParam->GetPadPitchWidth() +1.);
2183 PadNumber=(Int_t) (TMath::Sign((Float_t)PadNumber, y)+CentralPad);
2185 else continue; // electron out of sense wire range, lost at the sector edge
2187 Float_t aval = (absy<range-0.5) ? v(idx+4):v(idx+4)*IneffFactor;
2189 Float_t dist = y - (Float_t)(PadNumber-CentralPad)*fTPCParam->GetPadPitchWidth();
2191 for (Int_t i=0;i<7;i++){
2192 PadSignal[i]=fPRF2D->GetPRF(-dist+(i-3)*fTPCParam->GetPadPitchWidth(),xwire)*aval;
2194 PadSignal[i] *= fTPCParam->GetPadCoupling();
2198 Int_t LeftPad = TMath::Max(0,PadNumber-3);
2199 Int_t RightPad = TMath::Min(nPadsSignal-1,PadNumber+3);
2201 Int_t pmin=LeftPad-PadNumber+3; // lower index of the pad_signal vector
2202 Int_t pmax=RightPad-PadNumber+3; // upper index
2205 Float_t z_drift = z*zwidthm1;
2206 Float_t z_offset = z_drift-(Int_t)z_drift;
2208 Int_t FirstBucket = (Int_t)z_drift; // numbering starts from "0"
2210 for (Int_t i2=0;i2<4;i2++){
2211 Int_t TrueTime = FirstBucket+i2; // current time bucket
2212 Float_t dz = (Float_t(i2)+1.- z_offset)*zwidth;
2213 Float_t ampl = fRF->GetRF(dz);
2214 if((TrueTime>MAXTBKT-1)) break; // beyond the time range
2217 // loop over pads, from pmin to pmax
2219 for(Int_t i3=pmin;i3<pmax+1;i3++){
2220 Int_t TruePad = LeftPad+i3-pmin;
2222 if(TruePad<nPadsDiff || TruePad > nPadsSignal-nPadsDiff-1) continue;
2224 TruePad -= nPadsDiff;
2225 signal(TruePad,TrueTime)+=(PadSignal[i3]*ampl); // not converted to charge!
2227 } // end of loop over pads
2228 } // end of loop over time bins
2230 } // end of loop over electrons
2232 } // end of loop over tracks
2234 } // end of GetCrossTalk
2238 //_____________________________________________________________________________
2242 // Initialise TPC detector after definition of geometry
2247 for(i=0;i<35;i++) printf("*");
2248 printf(" TPC_INIT ");
2249 for(i=0;i<35;i++) printf("*");
2252 for(i=0;i<80;i++) printf("*");
2256 //_____________________________________________________________________________
2257 void AliTPC::MakeBranch(Option_t* option)
2260 // Create Tree branches for the TPC.
2262 Int_t buffersize = 4000;
2263 char branchname[10];
2264 sprintf(branchname,"%s",GetName());
2266 AliDetector::MakeBranch(option);
2268 char *D = strstr(option,"D");
2270 if (fDigits && gAlice->TreeD() && D) {
2271 gAlice->TreeD()->Branch(branchname,&fDigits, buffersize);
2272 printf("Making Branch %s for digits\n",branchname);
2275 char *R = strstr(option,"R");
2277 if (fClusters && gAlice->TreeR() && R) {
2278 gAlice->TreeR()->Branch(branchname,&fClusters, buffersize);
2279 printf("Making Branch %s for Clusters\n",branchname);
2283 //_____________________________________________________________________________
2284 void AliTPC::ResetDigits()
2287 // Reset number of digits and the digits array for this detector
2291 // if (fDigits) fDigits->Clear();
2292 if (fDigParam->GetArray()!=0) fDigParam->GetArray()->Clear();
2294 if (fClusters) fClusters->Clear();
2297 //_____________________________________________________________________________
2298 void AliTPC::SetSecAL(Int_t sec)
2300 //---------------------------------------------------
2301 // Activate/deactivate selection for lower sectors
2302 //---------------------------------------------------
2304 //-----------------------------------------------------------------
2305 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2306 //-----------------------------------------------------------------
2311 //_____________________________________________________________________________
2312 void AliTPC::SetSecAU(Int_t sec)
2314 //----------------------------------------------------
2315 // Activate/deactivate selection for upper sectors
2316 //---------------------------------------------------
2318 //-----------------------------------------------------------------
2319 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2320 //-----------------------------------------------------------------
2325 //_____________________________________________________________________________
2326 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2328 //----------------------------------------
2329 // Select active lower sectors
2330 //----------------------------------------
2332 //-----------------------------------------------------------------
2333 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2334 //-----------------------------------------------------------------
2344 //_____________________________________________________________________________
2345 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2346 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2347 Int_t s11 , Int_t s12)
2349 //--------------------------------
2350 // Select active upper sectors
2351 //--------------------------------
2353 //-----------------------------------------------------------------
2354 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2355 //-----------------------------------------------------------------
2371 //_____________________________________________________________________________
2372 void AliTPC::SetSens(Int_t sens)
2375 //-------------------------------------------------------------
2376 // Activates/deactivates the sensitive strips at the center of
2377 // the pad row -- this is for the space-point resolution calculations
2378 //-------------------------------------------------------------
2380 //-----------------------------------------------------------------
2381 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2382 //-----------------------------------------------------------------
2387 void AliTPC::SetSide(Float_t side)
2392 //____________________________________________________________________________
2393 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2394 Float_t p2,Float_t p3)
2409 //_____________________________________________________________________________
2410 void AliTPC::Streamer(TBuffer &R__b)
2413 // Stream an object of class AliTPC.
2415 if (R__b.IsReading()) {
2416 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
2417 AliDetector::Streamer(R__b);
2418 if (R__v < 2) return;
2422 fClustersIndex = new Int_t[fNsectors+1];
2423 fDigitsIndex = new Int_t[fNsectors+1];
2425 R__b.WriteVersion(AliTPC::IsA());
2426 AliDetector::Streamer(R__b);
2433 ClassImp(AliTPCcluster)
2435 //_____________________________________________________________________________
2436 AliTPCcluster::AliTPCcluster(Float_t *hits, Int_t *lab)
2439 // Creates a simulated cluster for the TPC
2441 fTracks[0] = lab[0];
2442 fTracks[1] = lab[1];
2443 fTracks[2] = lab[2];
2453 //_____________________________________________________________________________
2454 void AliTPCcluster::GetXYZ(Float_t *x, const AliTPCParam *par) const
2457 // Transformation from local to global coordinate system
2459 x[0]=par->GetPadRowRadii(fSector,fPadRow);
2462 par->CRXYZtoXYZ(x,fSector,fPadRow,1);
2466 //_____________________________________________________________________________
2467 Int_t AliTPCcluster::Compare(TObject * o)
2470 // compare two clusters according y coordinata
2472 AliTPCcluster *cl= (AliTPCcluster *)o;
2473 if (fY<cl->fY) return -1;
2474 if (fY==cl->fY) return 0;
2478 Bool_t AliTPCcluster::IsSortable() const
2481 //make AliTPCcluster sortabale
2488 ClassImp(AliTPCdigit)
2490 //_____________________________________________________________________________
2491 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2495 // Creates a TPC digit object
2497 fSector = digits[0];
2498 fPadRow = digits[1];
2501 fSignal = digits[4];
2507 //_____________________________________________________________________________
2508 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2512 // Creates a TPC hit object
2523 ClassImp(AliTPCtrack)
2525 //_____________________________________________________________________________
2526 AliTPCtrack::AliTPCtrack(Float_t *hits)
2529 // Default creator for a TPC reconstructed track object
2531 fX=hits[0]; // This is dummy code !
2534 AliTPCtrack::AliTPCtrack(const AliTPCcluster *c,const TVector& xx,
2535 const TMatrix& CC, Double_t xref, Double_t alpha):
2536 x(xx),C(CC),fClusters(200)
2538 //-----------------------------------------------------------------
2539 // This is the main track constructor.
2541 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2542 //-----------------------------------------------------------------
2546 fClusters.AddLast((AliTPCcluster*)(c));
2549 //_____________________________________________________________________________
2550 AliTPCtrack::AliTPCtrack(const AliTPCtrack& t) : x(t.x), C(t.C),
2551 fClusters(t.fClusters.GetEntriesFast())
2553 //-----------------------------------------------------------------
2554 // This is a track copy constructor.
2556 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2557 //-----------------------------------------------------------------
2561 int n=t.fClusters.GetEntriesFast();
2562 for (int i=0; i<n; i++) fClusters.AddLast(t.fClusters.UncheckedAt(i));
2565 //_____________________________________________________________________________
2566 Int_t AliTPCtrack::Compare(TObject *o) {
2567 //-----------------------------------------------------------------
2568 // This function compares tracks according to their curvature.
2570 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2571 //-----------------------------------------------------------------
2572 AliTPCtrack *t=(AliTPCtrack*)o;
2573 Double_t co=t->GetSigmaY2();
2574 Double_t c =GetSigmaY2();
2576 else if (c<co) return -1;
2580 //_____________________________________________________________________________
2581 int AliTPCtrack::PropagateTo(Double_t xk,Double_t x0,Double_t rho,Double_t pm)
2583 //-----------------------------------------------------------------
2584 // This function propagates a track to a reference plane x=xk.
2586 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2587 //-----------------------------------------------------------------
2588 if (TMath::Abs(x(2)*xk - x(3)) >= 0.999) {
2589 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Propagation failed !\n";
2593 Double_t x1=fX, x2=x1+0.5*(xk-x1), dx=x2-x1, y1=x(0), z1=x(1);
2594 Double_t c1=x(2)*x1 - x(3), r1=sqrt(1.- c1*c1);
2595 Double_t c2=x(2)*x2 - x(3), r2=sqrt(1.- c2*c2);
2597 x(0) += dx*(c1+c2)/(r1+r2);
2598 x(1) += dx*(c1+c2)/(c1*r2 + c2*r1)*x(4);
2600 TMatrix F(5,5); F.UnitMatrix();
2601 Double_t rr=r1+r2, cc=c1+c2, xx=x1+x2;
2602 F(0,2)= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
2603 F(0,3)=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
2604 Double_t cr=c1*r2+c2*r1;
2605 F(1,2)= dx*x(4)*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
2606 F(1,3)=-dx*x(4)*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
2608 TMatrix tmp(F,TMatrix::kMult,C);
2609 C.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2613 //Multiple scattering******************
2614 Double_t ey=x(2)*fX - x(3);
2615 Double_t ex=sqrt(1-ey*ey);
2617 TMatrix Q(5,5); Q=0.;
2618 Q(2,2)=ez*ez+ey*ey; Q(2,3)=-ex*ey; Q(2,4)=-ex*ez;
2619 Q(3,2)=Q(2,3); Q(3,3)= ez*ez+ex*ex; Q(3,4)=-ey*ez;
2620 Q(4,2)=Q(2,4); Q(4,3)= Q(3,4); Q(4,4)=1.;
2623 F(2,2)=-x(2)*ex; F(2,3)=-x(2)*ey;
2624 F(3,2)=-ex*(x(2)*fX-ey); F(3,3)=-(1.+ x(2)*fX*ey - ey*ey);
2625 F(4,2)=-ez*ex; F(4,3)=-ez*ey; F(4,4)=1.;
2628 Q.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2630 Double_t p2=GetPt()*GetPt()*(1.+x(4)*x(4));
2631 Double_t beta2=p2/(p2 + pm*pm);
2632 Double_t d=sqrt((x1-fX)*(x1-fX)+(y1-x(0))*(y1-x(0))+(z1-x(1))*(z1-x(1)));
2634 Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho;
2638 //Energy losses************************
2639 Double_t dE=0.153e-3/beta2*(log(5940*beta2/(1-beta2)) - beta2)*d*rho;
2640 if (x1 < x2) dE=-dE;
2641 x(2)*=(1.- sqrt(p2+pm*pm)/p2*dE);
2642 //x(3)*=(1.- sqrt(p2+pm*pm)/p2*dE);
2644 x1=fX; x2=xk; y1=x(0); z1=x(1);
2645 c1=x(2)*x1 - x(3); r1=sqrt(1.- c1*c1);
2646 c2=x(2)*x2 - x(3); r2=sqrt(1.- c2*c2);
2648 x(0) += dx*(c1+c2)/(r1+r2);
2649 x(1) += dx*(c1+c2)/(c1*r2 + c2*r1)*x(4);
2652 rr=r1+r2; cc=c1+c2; xx=x1+x2;
2653 F(0,2)= dx*(rr*xx + cc*(c1*x1/r1+c2*x2/r2))/(rr*rr);
2654 F(0,3)=-dx*(2*rr + cc*(c1/r1 + c2/r2))/(rr*rr);
2656 F(1,2)= dx*x(4)*(cr*xx-cc*(r1*x2-c2*c1*x1/r1+r2*x1-c1*c2*x2/r2))/(cr*cr);
2657 F(1,3)=-dx*x(4)*(2*cr + cc*(c2*c1/r1-r1 + c1*c2/r2-r2))/(cr*cr);
2660 C.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2667 //_____________________________________________________________________________
2668 void AliTPCtrack::PropagateToVertex(Double_t x0,Double_t rho,Double_t pm)
2670 //-----------------------------------------------------------------
2671 // This function propagates tracks to the "vertex".
2673 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2674 //-----------------------------------------------------------------
2675 Double_t c=x(2)*fX - x(3);
2676 Double_t tgf=-x(3)/(x(2)*x(0) + sqrt(1-c*c));
2677 Double_t snf=tgf/sqrt(1.+ tgf*tgf);
2678 Double_t xv=(x(3)+snf)/x(2);
2679 PropagateTo(xv,x0,rho,pm);
2682 //_____________________________________________________________________________
2683 void AliTPCtrack::Update(const AliTPCcluster *c, Double_t chisq)
2685 //-----------------------------------------------------------------
2686 // This function associates a clusters with this track.
2688 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2689 //-----------------------------------------------------------------
2690 TMatrix H(2,5); H.UnitMatrix();
2691 TMatrix Ht(TMatrix::kTransposed,H);
2692 TVector m(2); m(0)=c->fY; m(1)=c->fZ;
2693 TMatrix V(2,2); V(0,0)=c->fSigmaY2; V(0,1)=0.; V(1,0)=0.; V(1,1)=c->fSigmaZ2;
2695 TMatrix tmp(H,TMatrix::kMult,C);
2696 TMatrix R(tmp,TMatrix::kMult,Ht); R+=V;
2698 Double_t det=(Double_t)R(0,0)*R(1,1) - (Double_t)R(0,1)*R(1,0);
2699 R(0,1)=R(0,0); R(0,0)=R(1,1); R(1,1)=R(0,1);
2700 R(1,0)*=-1; R(0,1)=R(1,0);
2705 TMatrix K(C,TMatrix::kMult,Ht); K*=R;
2708 x*=H; x-=m; x*=-1; x*=K; x+=savex;
2709 if (TMath::Abs(x(2)*fX-x(3)) >= 0.999) {
2710 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Filtering failed !\n";
2716 C.Mult(K,tmp); C-=saveC; C*=-1;
2718 fClusters.AddLast((AliTPCcluster*)c);
2722 //_____________________________________________________________________________
2723 int AliTPCtrack::Rotate(Double_t alpha)
2725 //-----------------------------------------------------------------
2726 // This function rotates this track.
2728 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2729 //-----------------------------------------------------------------
2732 Double_t x1=fX, y1=x(0);
2733 Double_t ca=cos(alpha), sa=sin(alpha);
2734 Double_t r1=x(2)*fX - x(3);
2737 x(0)=-x1*sa + y1*ca;
2738 x(3)=x(3)*ca + (x(2)*y1 + sqrt(1.- r1*r1))*sa;
2740 Double_t r2=x(2)*fX - x(3);
2741 if (TMath::Abs(r2) >= 0.999) {
2742 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Rotation failed !\n";
2746 Double_t y0=x(0) + sqrt(1.- r2*r2)/x(2);
2747 if ((x(0)-y0)*x(2) >= 0.) {
2748 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Rotation failed !!!\n";
2752 TMatrix F(5,5); F.UnitMatrix();
2755 F(3,2)=(y1 - r1*x1/sqrt(1.- r1*r1))*sa;
2756 F(3,3)= ca + sa*r1/sqrt(1.- r1*r1);
2757 TMatrix tmp(F,TMatrix::kMult,C);
2758 // Double_t dy2=C(0,0);
2759 C.Mult(tmp,TMatrix(TMatrix::kTransposed,F));
2760 // C(0,0)+=dy2*sa*sa*r1*r1/(1.- r1*r1);
2761 // C(1,1)+=dy2*sa*sa*x(4)*x(4)/(1.- r1*r1);
2766 //_____________________________________________________________________________
2767 void AliTPCtrack::UseClusters() const
2769 //-----------------------------------------------------------------
2770 // This function marks clusters associated with this track.
2772 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2773 //-----------------------------------------------------------------
2774 int num_of_clusters=fClusters.GetEntriesFast();
2775 for (int i=0; i<num_of_clusters; i++) {
2776 //if (i<=14) continue;
2777 AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(i);
2782 //_____________________________________________________________________________
2783 Double_t AliTPCtrack::GetPredictedChi2(const AliTPCcluster *c) const
2785 //-----------------------------------------------------------------
2786 // This function calculates a predicted chi2 increment.
2788 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2789 //-----------------------------------------------------------------
2790 TMatrix H(2,5); H.UnitMatrix();
2791 TVector m(2); m(0)=c->fY; m(1)=c->fZ;
2792 TMatrix V(2,2); V(0,0)=c->fSigmaY2; V(0,1)=0.; V(1,0)=0.; V(1,1)=c->fSigmaZ2;
2793 TVector res=x; res*=H; res-=m; //res*=-1;
2794 TMatrix tmp(H,TMatrix::kMult,C);
2795 TMatrix R(tmp,TMatrix::kMult,TMatrix(TMatrix::kTransposed,H)); R+=V;
2797 Double_t det=(Double_t)R(0,0)*R(1,1) - (Double_t)R(0,1)*R(1,0);
2798 if (TMath::Abs(det) < 1.e-10) {
2799 if (*this>4) cerr<<*this<<" AliTPCtrack warning: Singular matrix !\n";
2802 R(0,1)=R(0,0); R(0,0)=R(1,1); R(1,1)=R(0,1);
2803 R(1,0)*=-1; R(0,1)=R(1,0);
2813 //_____________________________________________________________________________
2814 struct S { int lab; int max; };
2815 int AliTPCtrack::GetLabel(int nrows) const
2817 //-----------------------------------------------------------------
2818 // This function returns the track label. If label<0, this track is fake.
2820 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2821 //-----------------------------------------------------------------
2822 int num_of_clusters=fClusters.GetEntriesFast();
2823 S *s=new S[num_of_clusters];
2825 for (i=0; i<num_of_clusters; i++) s[i].lab=s[i].max=0;
2828 for (i=0; i<num_of_clusters; i++) {
2829 AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(i);
2830 lab=TMath::Abs(c->fTracks[0]);
2832 for (j=0; j<num_of_clusters; j++)
2833 if (s[j].lab==lab || s[j].max==0) break;
2839 for (i=0; i<num_of_clusters; i++)
2840 if (s[i].max>max) {max=s[i].max; lab=s[i].lab;}
2844 for (i=0; i<num_of_clusters; i++) {
2845 AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(i);
2846 if (TMath::Abs(c->fTracks[1]) == lab ||
2847 TMath::Abs(c->fTracks[2]) == lab ) max++;
2850 if (1.-float(max)/num_of_clusters > 0.10) return -lab;
2852 int tail=int(0.08*nrows);
2853 if (num_of_clusters < tail) return lab;
2856 for (i=1; i<=tail; i++) {
2857 AliTPCcluster *c=(AliTPCcluster*)fClusters.UncheckedAt(num_of_clusters-i);
2858 if (lab == TMath::Abs(c->fTracks[0]) ||
2859 lab == TMath::Abs(c->fTracks[1]) ||
2860 lab == TMath::Abs(c->fTracks[2])) max++;
2862 if (max < int(0.5*tail)) return -lab;
2867 //_____________________________________________________________________________
2868 void AliTPCtrack::GetPxPyPz(Double_t& px, Double_t& py, Double_t& pz) const
2870 //-----------------------------------------------------------------
2871 // This function returns reconstructed track momentum in the global system.
2873 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2874 //-----------------------------------------------------------------
2875 Double_t pt=TMath::Abs(GetPt()); // GeV/c
2876 Double_t r=x(2)*fX-x(3);
2877 Double_t y0=x(0) + sqrt(1.- r*r)/x(2);
2878 px=-pt*(x(0)-y0)*x(2); //cos(phi);
2879 py=-pt*(x(3)-fX*x(2)); //sin(phi);
2881 Double_t tmp=px*TMath::Cos(fAlpha) - py*TMath::Sin(fAlpha);
2882 py=px*TMath::Sin(fAlpha) + py*TMath::Cos(fAlpha);
2886 //_____________________________________________________________________________
2887 Double_t AliTPCtrack::GetdEdX(Double_t low, Double_t up) const {
2888 //-----------------------------------------------------------------
2889 // This funtion calculates dE/dX within the "low" and "up" cuts.
2891 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2892 //-----------------------------------------------------------------
2893 int ncl=fClusters.GetEntriesFast();
2895 Double_t *q=new Double_t[ncl];
2897 for (i=0; i<ncl; i++) {
2898 AliTPCcluster *cl=(AliTPCcluster*)(fClusters.UncheckedAt(i));
2899 // if (cl->fdEdX > 3000) continue;
2900 if (cl->fdEdX <= 0) continue;
2908 for (i=0; i<n-1; i++) {
2909 if (q[i]<=q[i+1]) continue;
2910 Double_t tmp=q[i]; q[i]=q[i+1]; q[i+1]=tmp;
2915 int nl=int(low*n), nu=int(up *n);
2917 for (i=nl; i<=nu; i++) dedx += q[i];
2922 //_________________________________________________________________________
2924 // Classes for internal tracking use
2925 //_________________________________________________________________________
2926 void AliTPCRow::InsertCluster(const AliTPCcluster* c) {
2927 //-----------------------------------------------------------------------
2928 // Insert a cluster into this pad row in accordence with its y-coordinate
2930 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2931 //-----------------------------------------------------------------------
2932 if (num_of_clusters==MAX_CLUSTER_PER_ROW) {
2933 cerr<<"AliTPCRow::InsertCluster(): Too many clusters !\n"; return;
2935 if (num_of_clusters==0) {clusters[num_of_clusters++]=c; return;}
2937 memmove(clusters+i+1 ,clusters+i,(num_of_clusters-i)*sizeof(AliTPCcluster*));
2938 clusters[i]=c; num_of_clusters++;
2941 int AliTPCRow::Find(Double_t y) const {
2942 //-----------------------------------------------------------------------
2943 // Return the index of the nearest cluster
2945 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
2946 //-----------------------------------------------------------------------
2947 if (y <= clusters[0]->fY) return 0;
2948 if (y > clusters[num_of_clusters-1]->fY) return num_of_clusters;
2949 int b=0, e=num_of_clusters-1, m=(b+e)/2;
2950 for (; b<e; m=(b+e)/2) {
2951 if (y > clusters[m]->fY) b=m+1;