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 **************************************************************************/
18 Revision 1.24 2000/10/05 16:06:09 kowal2
19 Forward declarations. Changes due to a new class AliComplexCluster.
21 Revision 1.23 2000/10/02 21:28:18 fca
22 Removal of useless dependecies via forward declarations
24 Revision 1.22 2000/07/10 20:57:39 hristov
25 Update of TPC code and macros by M.Kowalski
27 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
28 Changes to obey the coding rules
30 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
31 Splitted from AliTPCtracking
33 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
34 Changed parameter settings
36 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
38 Defaults loaded automatically (hard-wired)
39 Optional parameters can be set via macro called in the constructor
41 Revision 1.19 2000/04/18 19:00:59 fca
42 Small bug fixes to TPC files
44 Revision 1.18 2000/04/17 09:37:33 kowal2
45 removed obsolete AliTPCDigitsDisplay.C
47 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
49 New, experimental data structure from M. Ivanov
50 New tracking algorithm
51 Different pad geometry for different sectors
52 Digitization rewritten
54 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
55 Not used anymore - removed
57 Revision 1.17 2000/01/19 17:17:30 fca
58 Introducing a list of lists of hits -- more hits allowed for detector now
60 Revision 1.16 1999/11/05 09:29:23 fca
61 Accept only signals > 0
63 Revision 1.15 1999/10/08 06:26:53 fca
64 Removed ClustersIndex - not used anymore
66 Revision 1.14 1999/09/29 09:24:33 fca
67 Introduction of the Copyright and cvs Log
71 ///////////////////////////////////////////////////////////////////////////////
73 // Time Projection Chamber //
74 // This class contains the basic functions for the Time Projection Chamber //
75 // detector. Functions specific to one particular geometry are //
76 // contained in the derived classes //
80 <img src="picts/AliTPCClass.gif">
85 ///////////////////////////////////////////////////////////////////////////////
93 #include <TGeometry.h>
96 #include <TObjectTable.h>
97 #include "TParticle.h"
101 #include <iostream.h>
107 #include "AliTPCParamSR.h"
108 #include "AliTPCPRF2D.h"
109 #include "AliTPCRF1D.h"
110 #include "AliDigits.h"
111 #include "AliSimDigits.h"
112 #include "AliTPCTrackHits.h"
113 #include "AliPoints.h"
114 #include "AliArrayBranch.h"
117 #include "AliTPCDigitsArray.h"
118 #include "AliComplexCluster.h"
119 #include "AliClusters.h"
120 #include "AliTPCClustersRow.h"
121 #include "AliTPCClustersArray.h"
123 #include "AliTPCcluster.h"
124 #include "AliTPCclusterer.h"
125 #include "AliTPCtracker.h"
127 #include <TInterpreter.h>
134 //_____________________________________________________________________________
138 // Default constructor
153 //_____________________________________________________________________________
154 AliTPC::AliTPC(const char *name, const char *title)
155 : AliDetector(name,title)
158 // Standard constructor
162 // Initialise arrays of hits and digits
163 fHits = new TClonesArray("AliTPChit", 176);
164 gAlice->AddHitList(fHits);
169 fTrackHits = new AliTPCTrackHits; //MI - 13.09.2000
170 fTrackHits->SetHitPrecision(0.002);
171 fTrackHits->SetStepPrecision(0.003);
172 fTrackHits->SetMaxDistance(100);
176 // Initialise counters
182 // Initialise color attributes
183 SetMarkerColor(kYellow);
186 // Set TPC parameters
189 if (!strcmp(title,"Default")) {
190 fTPCParam = new AliTPCParamSR;
192 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
198 //_____________________________________________________________________________
209 delete fTrackHits; //MI 15.09.2000
212 //_____________________________________________________________________________
213 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
216 // Add a hit to the list
218 // TClonesArray &lhits = *fHits;
219 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
221 TClonesArray &lhits = *fHits;
222 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
225 AddHit2(track,vol,hits);
228 //_____________________________________________________________________________
229 void AliTPC::BuildGeometry()
233 // Build TPC ROOT TNode geometry for the event display
238 const int kColorTPC=19;
239 char name[5], title[25];
240 const Double_t kDegrad=TMath::Pi()/180;
241 const Double_t kRaddeg=180./TMath::Pi();
244 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
245 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
247 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
248 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
250 Int_t nLo = fTPCParam->GetNInnerSector()/2;
251 Int_t nHi = fTPCParam->GetNOuterSector()/2;
253 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
254 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
255 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
256 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
259 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
260 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
266 // Get ALICE top node
269 nTop=gAlice->GetGeometry()->GetNode("alice");
273 rl = fTPCParam->GetInnerRadiusLow();
274 ru = fTPCParam->GetInnerRadiusUp();
278 sprintf(name,"LS%2.2d",i);
280 sprintf(title,"TPC low sector %3d",i);
283 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
284 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
285 tubs->SetNumberOfDivisions(1);
287 nNode = new TNode(name,title,name,0,0,0,"");
288 nNode->SetLineColor(kColorTPC);
294 rl = fTPCParam->GetOuterRadiusLow();
295 ru = fTPCParam->GetOuterRadiusUp();
298 sprintf(name,"US%2.2d",i);
300 sprintf(title,"TPC upper sector %d",i);
302 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
303 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
304 tubs->SetNumberOfDivisions(1);
306 nNode = new TNode(name,title,name,0,0,0,"");
307 nNode->SetLineColor(kColorTPC);
313 //_____________________________________________________________________________
314 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
317 // Calculate distance from TPC to mouse on the display
323 void AliTPC::Clusters2Tracks(TFile *of) {
324 //-----------------------------------------------------------------
325 // This is a track finder.
326 //-----------------------------------------------------------------
327 AliTPCtracker::Clusters2Tracks(fTPCParam,of);
330 //_____________________________________________________________________________
331 void AliTPC::CreateMaterials()
333 //-----------------------------------------------
334 // Create Materials for for TPC simulations
335 //-----------------------------------------------
337 //-----------------------------------------------------------------
338 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
339 //-----------------------------------------------------------------
341 Int_t iSXFLD=gAlice->Field()->Integ();
342 Float_t sXMGMX=gAlice->Field()->Max();
344 Float_t amat[5]; // atomic numbers
345 Float_t zmat[5]; // z
346 Float_t wmat[5]; // proportions
352 //***************** Gases *************************
354 //-------------------------------------------------
356 //-------------------------------------------------
367 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
377 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
380 //--------------------------------------------------------------
382 //--------------------------------------------------------------
399 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
401 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
416 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
418 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
434 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
436 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
438 //----------------------------------------------------------------
439 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
440 //----------------------------------------------------------------
446 Float_t rho,absl,X0,buf[1];
450 for(nc = 0;nc<fNoComp;nc++)
453 // retrive material constants
455 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
460 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
462 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
463 density += fMixtProp[nc]*rho; // density of the mixture
467 // mixture proportions by weight!
469 for(nc = 0;nc<fNoComp;nc++)
472 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
474 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
475 apure[nnc] : amol[nnc])/am;
479 // Drift gases 1 - nonsensitive, 2 - sensitive
481 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
482 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
491 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
494 //----------------------------------------------------------------------
496 //----------------------------------------------------------------------
518 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
540 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
558 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
576 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
578 // G10 60% SiO2 + 40% epoxy, I use A and Z for SiO2
591 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz
593 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
595 AliMaterial(39,"G10",amat[0],zmat[0],density,999.,999.);
604 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
613 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
622 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
640 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
659 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
663 //----------------------------------------------------------
664 // tracking media for gases
665 //----------------------------------------------------------
667 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
668 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
669 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
670 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
672 //-----------------------------------------------------------
673 // tracking media for solids
674 //-----------------------------------------------------------
676 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
677 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
678 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
679 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
680 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
681 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
682 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
683 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
684 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
685 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
690 void AliTPC::Digits2Clusters(TFile *of)
692 //-----------------------------------------------------------------
693 // This is a simple cluster finder.
694 //-----------------------------------------------------------------
695 AliTPCclusterer::Digits2Clusters(fTPCParam,of);
698 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
699 extern Double_t SigmaZ2(Double_t, Double_t);
700 //_____________________________________________________________________________
701 void AliTPC::Hits2Clusters(TFile *of)
703 //--------------------------------------------------------
704 // TPC simple cluster generator from hits
705 // obtained from the TPC Fast Simulator
706 // The point errors are taken from the parametrization
707 //--------------------------------------------------------
709 //-----------------------------------------------------------------
710 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
711 //-----------------------------------------------------------------
712 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
713 // Jouri.Belikov@cern.ch
714 //----------------------------------------------------------------
716 /////////////////////////////////////////////////////////////////////////////
718 //---------------------------------------------------------------------
719 // ALICE TPC Cluster Parameters
720 //--------------------------------------------------------------------
724 // Cluster width in rphi
725 const Float_t kACrphi=0.18322;
726 const Float_t kBCrphi=0.59551e-3;
727 const Float_t kCCrphi=0.60952e-1;
728 // Cluster width in z
729 const Float_t kACz=0.19081;
730 const Float_t kBCz=0.55938e-3;
731 const Float_t kCCz=0.30428;
733 TDirectory *savedir=gDirectory;
736 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
741 printf("AliTPCParam MUST be created firstly\n");
745 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
747 TParticle *particle; // pointer to a given particle
748 AliTPChit *tpcHit; // pointer to a sigle TPC hit
749 TClonesArray *particles; //pointer to the particle list
753 Float_t pl,pt,tanth,rpad,ratio;
756 //---------------------------------------------------------------
757 // Get the access to the tracks
758 //---------------------------------------------------------------
760 TTree *tH = gAlice->TreeH();
761 Stat_t ntracks = tH->GetEntries();
762 particles=gAlice->Particles();
764 //Switch to the output file
767 fTPCParam->Write(fTPCParam->GetTitle());
768 AliTPCClustersArray carray;
769 carray.Setup(fTPCParam);
770 carray.SetClusterType("AliTPCcluster");
773 Int_t nclusters=0; //cluster counter
775 //------------------------------------------------------------
776 // Loop over all sectors (72 sectors for 20 deg
777 // segmentation for both lower and upper sectors)
778 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
779 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
781 // First cluster for sector 0 starts at "0"
782 //------------------------------------------------------------
784 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
786 fTPCParam->AdjustCosSin(isec,cph,sph);
788 //------------------------------------------------------------
790 //------------------------------------------------------------
792 for(Int_t track=0;track<ntracks;track++){
796 // Get number of the TPC hits
798 // nhits=fHits->GetEntriesFast();
801 tpcHit = (AliTPChit*)FirstHit(-1);
805 // for(Int_t hit=0;hit<nhits;hit++){
806 //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
810 if (tpcHit->fQ == 0.) {
811 tpcHit = (AliTPChit*) NextHit();
812 continue; //information about track (I.Belikov)
814 sector=tpcHit->fSector; // sector number
817 // if(sector != isec) continue; //terminate iteration
820 tpcHit = (AliTPChit*) NextHit();
823 ipart=tpcHit->Track();
824 particle=(TParticle*)particles->UncheckedAt(ipart);
827 if(pt < 1.e-9) pt=1.e-9;
829 tanth = TMath::Abs(tanth);
830 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
831 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
833 // space-point resolutions
835 sigmaRphi=SigmaY2(rpad,tanth,pt);
836 sigmaZ =SigmaZ2(rpad,tanth );
840 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
841 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
843 // temporary protection
845 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
846 if(sigmaZ < 0.) sigmaZ=0.4e-3;
847 if(clRphi < 0.) clRphi=2.5e-3;
848 if(clZ < 0.) clZ=2.5e-5;
853 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
854 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
856 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
857 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
858 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
859 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
860 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
861 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
862 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
863 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
864 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
865 xyz[2]=sigmaRphi; // fSigmaY2
866 xyz[3]=sigmaZ; // fSigmaZ2
867 xyz[4]=tpcHit->fQ; // q
869 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
870 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
872 Int_t tracks[3]={tpcHit->Track(), -1, -1};
873 AliTPCcluster cluster(tracks,xyz);
875 clrow->InsertCluster(&cluster); nclusters++;
877 tpcHit = (AliTPChit*)NextHit();
880 } // end of loop over hits
882 } // end of loop over tracks
884 Int_t nrows=fTPCParam->GetNRow(isec);
885 for (Int_t irow=0; irow<nrows; irow++) {
886 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
887 if (!clrow) continue;
888 carray.StoreRow(isec,irow);
889 carray.ClearRow(isec,irow);
892 } // end of loop over sectors
894 cerr<<"Number of made clusters : "<<nclusters<<" \n";
896 carray.GetTree()->Write();
898 savedir->cd(); //switch back to the input file
902 //_________________________________________________________________
903 void AliTPC::Hits2ExactClustersSector(Int_t isec)
905 //--------------------------------------------------------
906 //calculate exact cross point of track and given pad row
907 //resulting values are expressed in "digit" coordinata
908 //--------------------------------------------------------
910 //-----------------------------------------------------------------
911 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
912 //-----------------------------------------------------------------
914 if (fClustersArray==0){
918 TParticle *particle; // pointer to a given particle
919 AliTPChit *tpcHit; // pointer to a sigle TPC hit
920 TClonesArray *particles; //pointer to the particle list
923 const Int_t kcmaxhits=30000;
924 TVector * xxxx = new TVector(kcmaxhits*4);
925 TVector & xxx = *xxxx;
926 Int_t maxhits = kcmaxhits;
927 //construct array for each padrow
928 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
929 fClustersArray->CreateRow(isec,i);
931 //---------------------------------------------------------------
932 // Get the access to the tracks
933 //---------------------------------------------------------------
935 TTree *tH = gAlice->TreeH();
936 Stat_t ntracks = tH->GetEntries();
937 particles=gAlice->Particles();
938 Int_t npart = particles->GetEntriesFast();
940 //------------------------------------------------------------
942 //------------------------------------------------------------
944 for(Int_t track=0;track<ntracks;track++){
948 // Get number of the TPC hits and a pointer
951 nhits=fHits->GetEntriesFast();
955 Int_t currentIndex=0;
956 Int_t lastrow=-1; //last writen row
957 for(Int_t hit=0;hit<nhits;hit++){
958 tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
959 if (tpcHit==0) continue;
960 sector=tpcHit->fSector; // sector number
961 if(sector != isec) continue;
962 ipart=tpcHit->Track();
963 if (ipart<npart) particle=(TParticle*)particles->UncheckedAt(ipart);
967 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
968 Int_t index[3]={1,isec,0};
969 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
970 if (currentrow<0) continue;
971 if (lastrow<0) lastrow=currentrow;
972 if (currentrow==lastrow){
973 if ( currentIndex>=maxhits){
975 xxx.ResizeTo(4*maxhits);
977 xxx(currentIndex*4)=x[0];
978 xxx(currentIndex*4+1)=x[1];
979 xxx(currentIndex*4+2)=x[2];
980 xxx(currentIndex*4+3)=tpcHit->fQ;
996 for (Int_t index=0;index<currentIndex;index++){
998 x=x2=x3=x4=xxx(index*4);
1006 sumy+=xxx(index*4+1);
1007 sumxy+=xxx(index*4+1)*x;
1008 sumx2y+=xxx(index*4+1)*x2;
1009 sumz+=xxx(index*4+2);
1010 sumxz+=xxx(index*4+2)*x;
1011 sumx2z+=xxx(index*4+2)*x2;
1012 sumq+=xxx(index*4+3);
1014 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1015 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1016 sumx2*(sumx*sumx3-sumx2*sumx2);
1018 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1019 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1020 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1021 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1023 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1024 sumx2*(sumx*sumx2y-sumx2*sumxy);
1025 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1026 sumx2*(sumx*sumx2z-sumx2*sumxz);
1028 Float_t y=detay/det+centralPad;
1029 Float_t z=detaz/det;
1030 Float_t by=detby/det; //y angle
1031 Float_t bz=detbz/det; //z angle
1032 sumy/=Float_t(currentIndex);
1033 sumz/=Float_t(currentIndex);
1034 AliComplexCluster cl;
1040 cl.fTracks[0]=ipart;
1042 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1043 if (row!=0) row->InsertCluster(&cl);
1046 } //end of calculating cluster for given row
1050 } // end of loop over hits
1051 } // end of loop over tracks
1052 //write padrows to tree
1053 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1054 fClustersArray->StoreRow(isec,ii);
1055 fClustersArray->ClearRow(isec,ii);
1061 //__________________________________________________________________
1062 void AliTPC::Hits2Digits()
1064 //----------------------------------------------------
1065 // Loop over all sectors
1066 //----------------------------------------------------
1069 printf("AliTPCParam MUST be created firstly\n");
1073 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1078 //_____________________________________________________________________________
1079 void AliTPC::Hits2DigitsSector(Int_t isec)
1081 //-------------------------------------------------------------------
1082 // TPC conversion from hits to digits.
1083 //-------------------------------------------------------------------
1085 //-----------------------------------------------------------------
1086 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1087 //-----------------------------------------------------------------
1089 //-------------------------------------------------------
1090 // Get the access to the track hits
1091 //-------------------------------------------------------
1094 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1095 Stat_t ntracks = tH->GetEntries();
1099 //-------------------------------------------
1100 // Only if there are any tracks...
1101 //-------------------------------------------
1105 printf("*** Processing sector number %d ***\n",isec);
1107 Int_t nrows =fTPCParam->GetNRow(isec);
1109 row= new TObjArray* [nrows];
1111 MakeSector(isec,nrows,tH,ntracks,row);
1113 //--------------------------------------------------------
1114 // Digitize this sector, row by row
1115 // row[i] is the pointer to the TObjArray of TVectors,
1116 // each one containing electrons accepted on this
1117 // row, assigned into tracks
1118 //--------------------------------------------------------
1122 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree();
1124 for (i=0;i<nrows;i++){
1126 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1128 DigitizeRow(i,isec,row);
1130 fDigitsArray->StoreRow(isec,i);
1132 Int_t ndig = dig->GetDigitSize();
1134 printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1136 fDigitsArray->ClearRow(isec,i);
1139 } // end of the sector digitization
1141 for(i=0;i<nrows;i++){
1146 delete [] row; // delete the array of pointers to TObjArray-s
1150 } // end of Hits2DigitsSector
1153 //_____________________________________________________________________________
1154 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1156 //-----------------------------------------------------------
1157 // Single row digitization, coupling from the neighbouring
1158 // rows taken into account
1159 //-----------------------------------------------------------
1161 //-----------------------------------------------------------------
1162 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1163 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1164 //-----------------------------------------------------------------
1167 Float_t zerosup = fTPCParam->GetZeroSup();
1168 Int_t nrows =fTPCParam->GetNRow(isec);
1169 fCurrentIndex[1]= isec;
1172 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1173 Int_t nofTbins = fTPCParam->GetMaxTBin();
1174 Int_t indexRange[4];
1176 // Integrated signal for this row
1177 // and a single track signal
1179 TMatrix *m1 = new TMatrix(0,nofPads,0,nofTbins); // integrated
1180 TMatrix *m2 = new TMatrix(0,nofPads,0,nofTbins); // single
1182 TMatrix &total = *m1;
1184 // Array of pointers to the label-signal list
1186 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1187 Float_t **pList = new Float_t* [nofDigits];
1191 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1195 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1196 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1197 for (Int_t row= row1;row<=row2;row++){
1198 Int_t nTracks= rows[row]->GetEntries();
1199 for (i1=0;i1<nTracks;i1++){
1200 fCurrentIndex[2]= row;
1201 fCurrentIndex[3]=irow;
1203 m2->Zero(); // clear single track signal matrix
1204 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1205 GetList(trackLabel,nofPads,m2,indexRange,pList);
1207 else GetSignal(rows[row],i1,0,m1,indexRange);
1213 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1214 for(Int_t ip=0;ip<nofPads;ip++){
1215 for(Int_t it=0;it<nofTbins;it++){
1217 Float_t q = total(ip,it);
1219 Int_t gi =it*nofPads+ip; // global index
1221 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
1225 if(q <=zerosup) continue; // do not fill zeros
1226 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1229 // "real" signal or electronic noise (list = -1)?
1232 for(Int_t j1=0;j1<3;j1++){
1233 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1238 <A NAME="AliDigits"></A>
1239 using of AliDigits object
1242 dig->SetDigitFast((Short_t)q,it,ip);
1243 if (fDigitsArray->IsSimulated())
1245 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1246 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1247 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1251 } // end of loop over time buckets
1252 } // end of lop over pads
1255 // This row has been digitized, delete nonused stuff
1258 for(lp=0;lp<nofDigits;lp++){
1259 if(pList[lp]) delete [] pList[lp];
1268 } // end of DigitizeRow
1270 //_____________________________________________________________________________
1272 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, TMatrix *m1, TMatrix *m2,
1276 //---------------------------------------------------------------
1277 // Calculates 2-D signal (pad,time) for a single track,
1278 // returns a pointer to the signal matrix and the track label
1279 // No digitization is performed at this level!!!
1280 //---------------------------------------------------------------
1282 //-----------------------------------------------------------------
1283 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1284 // Modified: Marian Ivanov
1285 //-----------------------------------------------------------------
1289 tv = (TVector*)p1->At(ntr); // pointer to a track
1292 Float_t label = v(0);
1293 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1295 Int_t nElectrons = (tv->GetNrows()-1)/4;
1296 indexRange[0]=9999; // min pad
1297 indexRange[1]=-1; // max pad
1298 indexRange[2]=9999; //min time
1299 indexRange[3]=-1; // max time
1301 // Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1303 TMatrix &signal = *m1;
1304 TMatrix &total = *m2;
1306 // Loop over all electrons
1308 for(Int_t nel=0; nel<nElectrons; nel++){
1310 Float_t aval = v(idx+4);
1311 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1312 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1313 Int_t n = fTPCParam->CalcResponse(xyz,fCurrentIndex,fCurrentIndex[3]);
1315 if (n>0) for (Int_t i =0; i<n; i++){
1316 Int_t *index = fTPCParam->GetResBin(i);
1317 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1318 if ( ( pad<(fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]))) && (pad>0)) {
1319 Int_t time=index[2];
1320 Float_t weight = fTPCParam->GetResWeight(i); //we normalise response to ADC channel
1321 weight *= eltoadcfac;
1323 if (m1!=0) signal(pad,time)+=weight;
1324 total(pad,time)+=weight;
1325 indexRange[0]=TMath::Min(indexRange[0],pad);
1326 indexRange[1]=TMath::Max(indexRange[1],pad);
1327 indexRange[2]=TMath::Min(indexRange[2],time);
1328 indexRange[3]=TMath::Max(indexRange[3],time);
1331 } // end of loop over electrons
1333 return label; // returns track label when finished
1336 //_____________________________________________________________________________
1337 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *indexRange,
1340 //----------------------------------------------------------------------
1341 // Updates the list of tracks contributing to digits for a given row
1342 //----------------------------------------------------------------------
1344 //-----------------------------------------------------------------
1345 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1346 //-----------------------------------------------------------------
1348 TMatrix &signal = *m;
1350 // lop over nonzero digits
1352 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1353 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1356 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1358 if(signal(ip,it)<0.5) continue;
1361 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1363 if(!pList[globalIndex]){
1366 // Create new list (6 elements - 3 signals and 3 labels),
1369 pList[globalIndex] = new Float_t [6];
1373 *pList[globalIndex] = -1.;
1374 *(pList[globalIndex]+1) = -1.;
1375 *(pList[globalIndex]+2) = -1.;
1376 *(pList[globalIndex]+3) = -1.;
1377 *(pList[globalIndex]+4) = -1.;
1378 *(pList[globalIndex]+5) = -1.;
1381 *pList[globalIndex] = label;
1382 *(pList[globalIndex]+3) = signal(ip,it);
1386 // check the signal magnitude
1388 Float_t highest = *(pList[globalIndex]+3);
1389 Float_t middle = *(pList[globalIndex]+4);
1390 Float_t lowest = *(pList[globalIndex]+5);
1393 // compare the new signal with already existing list
1396 if(signal(ip,it)<lowest) continue; // neglect this track
1400 if (signal(ip,it)>highest){
1401 *(pList[globalIndex]+5) = middle;
1402 *(pList[globalIndex]+4) = highest;
1403 *(pList[globalIndex]+3) = signal(ip,it);
1405 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1406 *(pList[globalIndex]+1) = *pList[globalIndex];
1407 *pList[globalIndex] = label;
1409 else if (signal(ip,it)>middle){
1410 *(pList[globalIndex]+5) = middle;
1411 *(pList[globalIndex]+4) = signal(ip,it);
1413 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1414 *(pList[globalIndex]+1) = label;
1417 *(pList[globalIndex]+5) = signal(ip,it);
1418 *(pList[globalIndex]+2) = label;
1422 } // end of loop over pads
1423 } // end of loop over time bins
1428 //___________________________________________________________________
1429 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1430 Stat_t ntracks,TObjArray **row)
1433 //-----------------------------------------------------------------
1434 // Prepares the sector digitization, creates the vectors of
1435 // tracks for each row of this sector. The track vector
1436 // contains the track label and the position of electrons.
1437 //-----------------------------------------------------------------
1439 //-----------------------------------------------------------------
1440 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1441 //-----------------------------------------------------------------
1443 Float_t gasgain = fTPCParam->GetGasGain();
1447 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1450 if (fHitType&2) branch = TH->GetBranch("TPC2");
1451 else branch = TH->GetBranch("TPC");
1454 //----------------------------------------------
1455 // Create TObjArray-s, one for each row,
1456 // each TObjArray will store the TVectors
1457 // of electrons, one TVector per each track.
1458 //----------------------------------------------
1460 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
1461 TVector **tracks = new TVector* [nrows]; //pointers to the track vectors
1462 for(i=0; i<nrows; i++){
1463 row[i] = new TObjArray;
1470 //--------------------------------------------------------------------
1471 // Loop over tracks, the "track" contains the full history
1472 //--------------------------------------------------------------------
1474 Int_t previousTrack,currentTrack;
1475 previousTrack = -1; // nothing to store so far!
1477 for(Int_t track=0;track<ntracks;track++){
1478 Bool_t isInSector=kTRUE;
1483 TBranch * br = TH->GetBranch("fTrackHitsInfo");
1484 br->GetEvent(track);
1485 AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
1486 for (UInt_t j=0;j<ar->GetSize();j++){
1487 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
1490 if (!isInSector) continue;
1492 branch->GetEntry(track); // get next track
1496 tpcHit = (AliTPChit*)FirstHit(-1);
1498 //--------------------------------------------------------------
1500 //--------------------------------------------------------------
1505 Int_t sector=tpcHit->fSector; // sector number
1506 // if(sector != isec) continue;
1508 tpcHit = (AliTPChit*) NextHit();
1512 currentTrack = tpcHit->Track(); // track number
1515 if(currentTrack != previousTrack){
1517 // store already filled fTrack
1519 for(i=0;i<nrows;i++){
1520 if(previousTrack != -1){
1521 if(nofElectrons[i]>0){
1522 TVector &v = *tracks[i];
1523 v(0) = previousTrack;
1524 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1525 row[i]->Add(tracks[i]);
1528 delete tracks[i]; // delete empty TVector
1534 tracks[i] = new TVector(481); // TVectors for the next fTrack
1536 } // end of loop over rows
1538 previousTrack=currentTrack; // update track label
1541 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1543 //---------------------------------------------------
1544 // Calculate the electron attachment probability
1545 //---------------------------------------------------
1548 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
1549 /fTPCParam->GetDriftV();
1551 Float_t attProb = fTPCParam->GetAttCoef()*
1552 fTPCParam->GetOxyCont()*time; // fraction!
1554 //-----------------------------------------------
1555 // Loop over electrons
1556 //-----------------------------------------------
1559 for(Int_t nel=0;nel<qI;nel++){
1560 // skip if electron lost due to the attachment
1561 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1565 xyz[3]= (Float_t) (-gasgain*TMath::Log(gRandom->Rndm()));
1568 TransportElectron(xyz,index); //MI change -august
1570 fTPCParam->GetPadRow(xyz,index); //MI change august
1571 rowNumber = index[2];
1572 //transform position to local digit coordinates
1573 //relative to nearest pad row
1574 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
1575 nofElectrons[rowNumber]++;
1576 //----------------------------------
1577 // Expand vector if necessary
1578 //----------------------------------
1579 if(nofElectrons[rowNumber]>120){
1580 Int_t range = tracks[rowNumber]->GetNrows();
1581 if((nofElectrons[rowNumber])>(range-1)/4){
1583 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
1587 TVector &v = *tracks[rowNumber];
1588 Int_t idx = 4*nofElectrons[rowNumber]-3;
1590 v(idx)= xyz[0]; // X - pad row coordinate
1591 v(idx+1)=xyz[1]; // Y - pad coordinate (along the pad-row)
1592 v(idx+2)=xyz[2]; // Z - time bin coordinate
1593 v(idx+3)=xyz[3]; // avalanche size
1594 } // end of loop over electrons
1596 tpcHit = (AliTPChit*)NextHit();
1598 } // end of loop over hits
1599 } // end of loop over tracks
1602 // store remaining track (the last one) if not empty
1605 for(i=0;i<nrows;i++){
1606 if(nofElectrons[i]>0){
1607 TVector &v = *tracks[i];
1608 v(0) = previousTrack;
1609 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1610 row[i]->Add(tracks[i]);
1619 delete [] nofElectrons;
1622 } // end of MakeSector
1625 //_____________________________________________________________________________
1629 // Initialise TPC detector after definition of geometry
1634 for(i=0;i<35;i++) printf("*");
1635 printf(" TPC_INIT ");
1636 for(i=0;i<35;i++) printf("*");
1639 for(i=0;i<80;i++) printf("*");
1643 //_____________________________________________________________________________
1644 void AliTPC::MakeBranch(Option_t* option)
1647 // Create Tree branches for the TPC.
1649 Int_t buffersize = 4000;
1650 char branchname[10];
1651 sprintf(branchname,"%s",GetName());
1653 AliDetector::MakeBranch(option);
1655 char *d = strstr(option,"D");
1657 if (fDigits && gAlice->TreeD() && d) {
1658 gAlice->TreeD()->Branch(branchname,&fDigits, buffersize);
1659 printf("Making Branch %s for digits\n",branchname);
1662 if (fHitType&2) MakeBranch2(option); // MI change 14.09.2000
1665 //_____________________________________________________________________________
1666 void AliTPC::ResetDigits()
1669 // Reset number of digits and the digits array for this detector
1672 if (fDigits) fDigits->Clear();
1675 //_____________________________________________________________________________
1676 void AliTPC::SetSecAL(Int_t sec)
1678 //---------------------------------------------------
1679 // Activate/deactivate selection for lower sectors
1680 //---------------------------------------------------
1682 //-----------------------------------------------------------------
1683 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1684 //-----------------------------------------------------------------
1689 //_____________________________________________________________________________
1690 void AliTPC::SetSecAU(Int_t sec)
1692 //----------------------------------------------------
1693 // Activate/deactivate selection for upper sectors
1694 //---------------------------------------------------
1696 //-----------------------------------------------------------------
1697 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1698 //-----------------------------------------------------------------
1703 //_____________________________________________________________________________
1704 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
1706 //----------------------------------------
1707 // Select active lower sectors
1708 //----------------------------------------
1710 //-----------------------------------------------------------------
1711 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1712 //-----------------------------------------------------------------
1722 //_____________________________________________________________________________
1723 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
1724 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
1725 Int_t s11 , Int_t s12)
1727 //--------------------------------
1728 // Select active upper sectors
1729 //--------------------------------
1731 //-----------------------------------------------------------------
1732 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1733 //-----------------------------------------------------------------
1749 //_____________________________________________________________________________
1750 void AliTPC::SetSens(Int_t sens)
1753 //-------------------------------------------------------------
1754 // Activates/deactivates the sensitive strips at the center of
1755 // the pad row -- this is for the space-point resolution calculations
1756 //-------------------------------------------------------------
1758 //-----------------------------------------------------------------
1759 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1760 //-----------------------------------------------------------------
1766 void AliTPC::SetSide(Float_t side=0.)
1768 // choice of the TPC side
1773 //____________________________________________________________________________
1774 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
1775 Float_t p2,Float_t p3)
1778 // gax mixture definition
1792 //_____________________________________________________________________________
1794 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
1797 // electron transport taking into account:
1799 // 2.ExB at the wires
1800 // 3. nonisochronity
1802 // xyz and index must be already transformed to system 1
1805 fTPCParam->Transform1to2(xyz,index);
1808 Float_t driftl=xyz[2];
1809 if(driftl<0.01) driftl=0.01;
1810 driftl=TMath::Sqrt(driftl);
1811 Float_t sigT = driftl*(fTPCParam->GetDiffT());
1812 Float_t sigL = driftl*(fTPCParam->GetDiffL());
1813 xyz[0]=gRandom->Gaus(xyz[0],sigT);
1814 xyz[1]=gRandom->Gaus(xyz[1],sigT);
1815 xyz[2]=gRandom->Gaus(xyz[2],sigL);
1819 if (fTPCParam->GetMWPCReadout()==kTRUE){
1821 fTPCParam->Transform2to2NearestWire(xyz,index);
1822 Float_t dx=xyz[0]-x1;
1823 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
1825 //add nonisochronity (not implemented yet)
1828 //_____________________________________________________________________________
1829 void AliTPC::Streamer(TBuffer &R__b)
1832 // Stream an object of class AliTPC.
1834 if (R__b.IsReading()) {
1835 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1836 AliDetector::Streamer(R__b);
1838 if (R__v < 2) return;
1843 R__b.WriteVersion(AliTPC::IsA());
1844 AliDetector::Streamer(R__b);
1851 ClassImp(AliTPCdigit)
1853 //_____________________________________________________________________________
1854 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
1858 // Creates a TPC digit object
1860 fSector = digits[0];
1861 fPadRow = digits[1];
1864 fSignal = digits[4];
1870 //_____________________________________________________________________________
1871 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
1875 // Creates a TPC hit object
1886 //________________________________________________________________________
1887 // Additional code because of the AliTPCTrackHits
1889 void AliTPC::MakeBranch2(Option_t *option)
1892 // Create a new branch in the current Root Tree
1893 // The branch of fHits is automatically split
1894 // MI change 14.09.2000
1895 if (fHitType&2==0) return;
1896 char branchname[10];
1897 sprintf(branchname,"%s2",GetName());
1899 // Get the pointer to the header
1900 char *cH = strstr(option,"H");
1902 if (fTrackHits && gAlice->TreeH() && cH) {
1903 // gAlice->TreeH()->Branch(branchname,&fTrackHits, fBufferSize);
1904 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHits,
1905 gAlice->TreeH(),fBufferSize,1);
1906 gAlice->TreeH()->GetListOfBranches()->Add(branch);
1907 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
1911 void AliTPC::SetTreeAddress()
1913 if (fHitType&1) AliDetector::SetTreeAddress();
1914 if (fHitType&2) SetTreeAddress2();
1917 void AliTPC::SetTreeAddress2()
1920 // Set branch address for the TrackHits Tree
1923 char branchname[20];
1924 sprintf(branchname,"%s2",GetName());
1926 // Branch address for hit tree
1927 TTree *treeH = gAlice->TreeH();
1929 branch = treeH->GetBranch(branchname);
1930 if (branch) branch->SetAddress(&fTrackHits);
1934 void AliTPC::FinishPrimary()
1936 if (fTrackHits) fTrackHits->FlushHitStack();
1940 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
1943 // add hit to the list
1944 TClonesArray &particles = *(gAlice->Particles());
1947 int primary = gAlice->GetPrimary(track);
1948 ((TParticle *)particles[primary])->SetBit(kKeepBit);
1952 gAlice->FlagTrack(track);
1954 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
1955 //if (hit->fTrack!=rtrack)
1956 // cout<<"bad track number\n";
1958 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
1959 hits[1],hits[2],(Int_t)hits[3]);
1962 void AliTPC::ResetHits()
1964 if (fHitType&1) AliDetector::ResetHits();
1965 if (fHitType&2) ResetHits2();
1968 void AliTPC::ResetHits2()
1972 if (fTrackHits) fTrackHits->Clear();
1975 AliHit* AliTPC::FirstHit(Int_t track)
1977 if (fHitType&2) return FirstHit2(track);
1978 return AliDetector::FirstHit(track);
1980 AliHit* AliTPC::NextHit()
1982 if (fHitType&2) return NextHit2();
1983 return AliDetector::NextHit();
1986 AliHit* AliTPC::FirstHit2(Int_t track)
1989 // Initialise the hit iterator
1990 // Return the address of the first hit for track
1991 // If track>=0 the track is read from disk
1992 // while if track<0 the first hit of the current
1993 // track is returned
1996 gAlice->ResetHits();
1997 gAlice->TreeH()->GetEvent(track);
2001 fTrackHits->First();
2002 return fTrackHits->GetHit();
2007 AliHit* AliTPC::NextHit2()
2010 //Return the next hit for the current track
2014 return fTrackHits->GetHit();
2020 void AliTPC::LoadPoints(Int_t)
2030 void AliTPC::RemapTrackHitIDs(Int_t *map)
2032 if (!fTrackHits) return;
2033 AliObjectArray * arr = fTrackHits->fTrackHitsInfo;
2034 for (UInt_t i=0;i<arr->GetSize();i++){
2035 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2036 info->fTrackID = map[info->fTrackID];
2042 //_____________________________________________________________________________
2043 void AliTPC::LoadPoints2(Int_t)
2046 // Store x, y, z of all hits in memory
2048 if (fTrackHits == 0) return;
2050 Int_t nhits = fTrackHits->GetEntriesFast();
2051 if (nhits == 0) return;
2052 Int_t tracks = gAlice->GetNtrack();
2053 if (fPoints == 0) fPoints = new TObjArray(tracks);
2056 Int_t *ntrk=new Int_t[tracks];
2057 Int_t *limi=new Int_t[tracks];
2058 Float_t **coor=new Float_t*[tracks];
2059 for(Int_t i=0;i<tracks;i++) {
2065 AliPoints *points = 0;
2068 Int_t chunk=nhits/4+1;
2070 // Loop over all the hits and store their position
2072 ahit = FirstHit2(-1);
2073 //for (Int_t hit=0;hit<nhits;hit++) {
2075 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2076 trk=ahit->GetTrack();
2077 if(ntrk[trk]==limi[trk]) {
2079 // Initialise a new track
2080 fp=new Float_t[3*(limi[trk]+chunk)];
2082 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2083 delete [] coor[trk];
2090 fp[3*ntrk[trk] ] = ahit->X();
2091 fp[3*ntrk[trk]+1] = ahit->Y();
2092 fp[3*ntrk[trk]+2] = ahit->Z();
2097 for(trk=0; trk<tracks; ++trk) {
2099 points = new AliPoints();
2100 points->SetMarkerColor(GetMarkerColor());
2101 points->SetMarkerSize(GetMarkerSize());
2102 points->SetDetector(this);
2103 points->SetParticle(trk);
2104 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2105 fPoints->AddAt(points,trk);
2106 delete [] coor[trk];
2116 //_____________________________________________________________________________
2117 void AliTPC::LoadPoints3(Int_t)
2120 // Store x, y, z of all hits in memory
2121 // - only intersection point with pad row
2122 if (fTrackHits == 0) return;
2124 Int_t nhits = fTrackHits->GetEntriesFast();
2125 if (nhits == 0) return;
2126 Int_t tracks = gAlice->GetNtrack();
2127 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2128 fPoints->Expand(2*tracks);
2131 Int_t *ntrk=new Int_t[tracks];
2132 Int_t *limi=new Int_t[tracks];
2133 Float_t **coor=new Float_t*[tracks];
2134 for(Int_t i=0;i<tracks;i++) {
2140 AliPoints *points = 0;
2143 Int_t chunk=nhits/4+1;
2145 // Loop over all the hits and store their position
2147 ahit = FirstHit2(-1);
2148 //for (Int_t hit=0;hit<nhits;hit++) {
2152 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2153 trk=ahit->GetTrack();
2154 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2155 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2156 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2157 if (currentrow!=lastrow){
2158 lastrow = currentrow;
2159 //later calculate intersection point
2160 if(ntrk[trk]==limi[trk]) {
2162 // Initialise a new track
2163 fp=new Float_t[3*(limi[trk]+chunk)];
2165 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2166 delete [] coor[trk];
2173 fp[3*ntrk[trk] ] = ahit->X();
2174 fp[3*ntrk[trk]+1] = ahit->Y();
2175 fp[3*ntrk[trk]+2] = ahit->Z();
2182 for(trk=0; trk<tracks; ++trk) {
2184 points = new AliPoints();
2185 points->SetMarkerColor(GetMarkerColor()+1);
2186 points->SetMarkerStyle(5);
2187 points->SetMarkerSize(0.2);
2188 points->SetDetector(this);
2189 points->SetParticle(trk);
2190 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2191 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2192 fPoints->AddAt(points,tracks+trk);
2193 delete [] coor[trk];
2204 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2208 //fill clones array with intersection of current point with the
2213 const Int_t kcmaxhits=30000;
2214 TVector * xxxx = new TVector(kcmaxhits*4);
2215 TVector & xxx = *xxxx;
2216 Int_t maxhits = kcmaxhits;
2219 AliTPChit * tpcHit=0;
2220 tpcHit = (AliTPChit*)FirstHit2(-1);
2221 Int_t currentIndex=0;
2222 Int_t lastrow=-1; //last writen row
2225 if (tpcHit==0) continue;
2226 sector=tpcHit->fSector; // sector number
2227 ipart=tpcHit->Track();
2231 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2232 Int_t index[3]={1,sector,0};
2233 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2234 if (currentrow<0) continue;
2235 if (lastrow<0) lastrow=currentrow;
2236 if (currentrow==lastrow){
2237 if ( currentIndex>=maxhits){
2239 xxx.ResizeTo(4*maxhits);
2241 xxx(currentIndex*4)=x[0];
2242 xxx(currentIndex*4+1)=x[1];
2243 xxx(currentIndex*4+2)=x[2];
2244 xxx(currentIndex*4+3)=tpcHit->fQ;
2248 if (currentIndex>2){
2260 for (Int_t index=0;index<currentIndex;index++){
2262 x=x2=x3=x4=xxx(index*4);
2270 sumy+=xxx(index*4+1);
2271 sumxy+=xxx(index*4+1)*x;
2272 sumx2y+=xxx(index*4+1)*x2;
2273 sumz+=xxx(index*4+2);
2274 sumxz+=xxx(index*4+2)*x;
2275 sumx2z+=xxx(index*4+2)*x2;
2276 sumq+=xxx(index*4+3);
2278 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2279 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2280 sumx2*(sumx*sumx3-sumx2*sumx2);
2282 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2283 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2284 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2285 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2287 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2288 sumx2*(sumx*sumx2y-sumx2*sumxy);
2289 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2290 sumx2*(sumx*sumx2z-sumx2*sumxz);
2292 Float_t y=detay/det+centralPad;
2293 Float_t z=detaz/det;
2294 Float_t by=detby/det; //y angle
2295 Float_t bz=detbz/det; //z angle
2296 sumy/=Float_t(currentIndex);
2297 sumz/=Float_t(currentIndex);
2299 AliComplexCluster cl;
2305 cl.fTracks[0]=ipart;
2307 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2308 if (row!=0) row->InsertCluster(&cl);
2311 } //end of calculating cluster for given row
2313 } // end of loop over hits