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.26 2001/01/10 07:59:43 kowal2
19 Corrections to load points from the noncompressed hits.
21 Revision 1.25 2000/11/02 07:25:31 kowal2
22 Changes due to the new hit structure.
25 Revision 1.24 2000/10/05 16:06:09 kowal2
26 Forward declarations. Changes due to a new class AliComplexCluster.
28 Revision 1.23 2000/10/02 21:28:18 fca
29 Removal of useless dependecies via forward declarations
31 Revision 1.22 2000/07/10 20:57:39 hristov
32 Update of TPC code and macros by M.Kowalski
34 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
35 Changes to obey the coding rules
37 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
38 Splitted from AliTPCtracking
40 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
41 Changed parameter settings
43 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
45 Defaults loaded automatically (hard-wired)
46 Optional parameters can be set via macro called in the constructor
48 Revision 1.19 2000/04/18 19:00:59 fca
49 Small bug fixes to TPC files
51 Revision 1.18 2000/04/17 09:37:33 kowal2
52 removed obsolete AliTPCDigitsDisplay.C
54 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
56 New, experimental data structure from M. Ivanov
57 New tracking algorithm
58 Different pad geometry for different sectors
59 Digitization rewritten
61 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
62 Not used anymore - removed
64 Revision 1.17 2000/01/19 17:17:30 fca
65 Introducing a list of lists of hits -- more hits allowed for detector now
67 Revision 1.16 1999/11/05 09:29:23 fca
68 Accept only signals > 0
70 Revision 1.15 1999/10/08 06:26:53 fca
71 Removed ClustersIndex - not used anymore
73 Revision 1.14 1999/09/29 09:24:33 fca
74 Introduction of the Copyright and cvs Log
78 ///////////////////////////////////////////////////////////////////////////////
80 // Time Projection Chamber //
81 // This class contains the basic functions for the Time Projection Chamber //
82 // detector. Functions specific to one particular geometry are //
83 // contained in the derived classes //
87 <img src="picts/AliTPCClass.gif">
92 ///////////////////////////////////////////////////////////////////////////////
100 #include <TGeometry.h>
103 #include <TObjectTable.h>
104 #include "TParticle.h"
108 #include <iostream.h>
114 #include "AliTPCParamSR.h"
115 #include "AliTPCPRF2D.h"
116 #include "AliTPCRF1D.h"
117 #include "AliDigits.h"
118 #include "AliSimDigits.h"
119 #include "AliTPCTrackHits.h"
120 #include "AliPoints.h"
121 #include "AliArrayBranch.h"
124 #include "AliTPCDigitsArray.h"
125 #include "AliComplexCluster.h"
126 #include "AliClusters.h"
127 #include "AliTPCClustersRow.h"
128 #include "AliTPCClustersArray.h"
130 #include "AliTPCcluster.h"
131 #include "AliTPCclusterer.h"
132 #include "AliTPCtracker.h"
134 #include <TInterpreter.h>
141 //_____________________________________________________________________________
145 // Default constructor
160 //_____________________________________________________________________________
161 AliTPC::AliTPC(const char *name, const char *title)
162 : AliDetector(name,title)
165 // Standard constructor
169 // Initialise arrays of hits and digits
170 fHits = new TClonesArray("AliTPChit", 176);
171 gAlice->AddHitList(fHits);
176 fTrackHits = new AliTPCTrackHits; //MI - 13.09.2000
177 fTrackHits->SetHitPrecision(0.002);
178 fTrackHits->SetStepPrecision(0.003);
179 fTrackHits->SetMaxDistance(100);
183 // Initialise counters
189 // Initialise color attributes
190 SetMarkerColor(kYellow);
193 // Set TPC parameters
196 if (!strcmp(title,"Default")) {
197 fTPCParam = new AliTPCParamSR;
199 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
205 //_____________________________________________________________________________
216 delete fTrackHits; //MI 15.09.2000
219 //_____________________________________________________________________________
220 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
223 // Add a hit to the list
225 // TClonesArray &lhits = *fHits;
226 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
228 TClonesArray &lhits = *fHits;
229 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
232 AddHit2(track,vol,hits);
235 //_____________________________________________________________________________
236 void AliTPC::BuildGeometry()
240 // Build TPC ROOT TNode geometry for the event display
245 const int kColorTPC=19;
246 char name[5], title[25];
247 const Double_t kDegrad=TMath::Pi()/180;
248 const Double_t kRaddeg=180./TMath::Pi();
251 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
252 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
254 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
255 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
257 Int_t nLo = fTPCParam->GetNInnerSector()/2;
258 Int_t nHi = fTPCParam->GetNOuterSector()/2;
260 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
261 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
262 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
263 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
266 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
267 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
273 // Get ALICE top node
276 nTop=gAlice->GetGeometry()->GetNode("alice");
280 rl = fTPCParam->GetInnerRadiusLow();
281 ru = fTPCParam->GetInnerRadiusUp();
285 sprintf(name,"LS%2.2d",i);
287 sprintf(title,"TPC low sector %3d",i);
290 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
291 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
292 tubs->SetNumberOfDivisions(1);
294 nNode = new TNode(name,title,name,0,0,0,"");
295 nNode->SetLineColor(kColorTPC);
301 rl = fTPCParam->GetOuterRadiusLow();
302 ru = fTPCParam->GetOuterRadiusUp();
305 sprintf(name,"US%2.2d",i);
307 sprintf(title,"TPC upper sector %d",i);
309 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
310 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
311 tubs->SetNumberOfDivisions(1);
313 nNode = new TNode(name,title,name,0,0,0,"");
314 nNode->SetLineColor(kColorTPC);
320 //_____________________________________________________________________________
321 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
324 // Calculate distance from TPC to mouse on the display
330 void AliTPC::Clusters2Tracks(TFile *of) {
331 //-----------------------------------------------------------------
332 // This is a track finder.
333 //-----------------------------------------------------------------
334 AliTPCtracker::Clusters2Tracks(fTPCParam,of);
337 //_____________________________________________________________________________
338 void AliTPC::CreateMaterials()
340 //-----------------------------------------------
341 // Create Materials for for TPC simulations
342 //-----------------------------------------------
344 //-----------------------------------------------------------------
345 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
346 //-----------------------------------------------------------------
348 Int_t iSXFLD=gAlice->Field()->Integ();
349 Float_t sXMGMX=gAlice->Field()->Max();
351 Float_t amat[5]; // atomic numbers
352 Float_t zmat[5]; // z
353 Float_t wmat[5]; // proportions
359 //***************** Gases *************************
361 //-------------------------------------------------
363 //-------------------------------------------------
374 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
384 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
387 //--------------------------------------------------------------
389 //--------------------------------------------------------------
406 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
408 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
423 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
425 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
441 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
443 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
445 //----------------------------------------------------------------
446 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
447 //----------------------------------------------------------------
453 Float_t rho,absl,X0,buf[1];
457 for(nc = 0;nc<fNoComp;nc++)
460 // retrive material constants
462 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
467 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
469 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
470 density += fMixtProp[nc]*rho; // density of the mixture
474 // mixture proportions by weight!
476 for(nc = 0;nc<fNoComp;nc++)
479 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
481 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
482 apure[nnc] : amol[nnc])/am;
486 // Drift gases 1 - nonsensitive, 2 - sensitive
488 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
489 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
498 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
501 //----------------------------------------------------------------------
503 //----------------------------------------------------------------------
525 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
547 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
565 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
583 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
585 // G10 60% SiO2 + 40% epoxy, I use A and Z for SiO2
598 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz
600 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
602 AliMaterial(39,"G10",amat[0],zmat[0],density,999.,999.);
611 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
620 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
629 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
647 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
666 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
670 //----------------------------------------------------------
671 // tracking media for gases
672 //----------------------------------------------------------
674 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
675 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
676 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
677 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
679 //-----------------------------------------------------------
680 // tracking media for solids
681 //-----------------------------------------------------------
683 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
684 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
685 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
686 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
687 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
688 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
689 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
690 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
691 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
692 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
697 void AliTPC::Digits2Clusters(TFile *of)
699 //-----------------------------------------------------------------
700 // This is a simple cluster finder.
701 //-----------------------------------------------------------------
702 AliTPCclusterer::Digits2Clusters(fTPCParam,of);
705 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
706 extern Double_t SigmaZ2(Double_t, Double_t);
707 //_____________________________________________________________________________
708 void AliTPC::Hits2Clusters(TFile *of)
710 //--------------------------------------------------------
711 // TPC simple cluster generator from hits
712 // obtained from the TPC Fast Simulator
713 // The point errors are taken from the parametrization
714 //--------------------------------------------------------
716 //-----------------------------------------------------------------
717 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
718 //-----------------------------------------------------------------
719 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
720 // Jouri.Belikov@cern.ch
721 //----------------------------------------------------------------
723 /////////////////////////////////////////////////////////////////////////////
725 //---------------------------------------------------------------------
726 // ALICE TPC Cluster Parameters
727 //--------------------------------------------------------------------
731 // Cluster width in rphi
732 const Float_t kACrphi=0.18322;
733 const Float_t kBCrphi=0.59551e-3;
734 const Float_t kCCrphi=0.60952e-1;
735 // Cluster width in z
736 const Float_t kACz=0.19081;
737 const Float_t kBCz=0.55938e-3;
738 const Float_t kCCz=0.30428;
740 TDirectory *savedir=gDirectory;
743 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
748 printf("AliTPCParam MUST be created firstly\n");
752 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
754 TParticle *particle; // pointer to a given particle
755 AliTPChit *tpcHit; // pointer to a sigle TPC hit
756 TClonesArray *particles; //pointer to the particle list
760 Float_t pl,pt,tanth,rpad,ratio;
763 //---------------------------------------------------------------
764 // Get the access to the tracks
765 //---------------------------------------------------------------
767 TTree *tH = gAlice->TreeH();
768 Stat_t ntracks = tH->GetEntries();
769 particles=gAlice->Particles();
771 //Switch to the output file
774 fTPCParam->Write(fTPCParam->GetTitle());
775 AliTPCClustersArray carray;
776 carray.Setup(fTPCParam);
777 carray.SetClusterType("AliTPCcluster");
780 Int_t nclusters=0; //cluster counter
782 //------------------------------------------------------------
783 // Loop over all sectors (72 sectors for 20 deg
784 // segmentation for both lower and upper sectors)
785 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
786 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
788 // First cluster for sector 0 starts at "0"
789 //------------------------------------------------------------
791 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
793 fTPCParam->AdjustCosSin(isec,cph,sph);
795 //------------------------------------------------------------
797 //------------------------------------------------------------
799 for(Int_t track=0;track<ntracks;track++){
803 // Get number of the TPC hits
805 // nhits=fHits->GetEntriesFast();
808 tpcHit = (AliTPChit*)FirstHit(-1);
812 // for(Int_t hit=0;hit<nhits;hit++){
813 //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
817 if (tpcHit->fQ == 0.) {
818 tpcHit = (AliTPChit*) NextHit();
819 continue; //information about track (I.Belikov)
821 sector=tpcHit->fSector; // sector number
824 // if(sector != isec) continue; //terminate iteration
827 tpcHit = (AliTPChit*) NextHit();
830 ipart=tpcHit->Track();
831 particle=(TParticle*)particles->UncheckedAt(ipart);
834 if(pt < 1.e-9) pt=1.e-9;
836 tanth = TMath::Abs(tanth);
837 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
838 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
840 // space-point resolutions
842 sigmaRphi=SigmaY2(rpad,tanth,pt);
843 sigmaZ =SigmaZ2(rpad,tanth );
847 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
848 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
850 // temporary protection
852 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
853 if(sigmaZ < 0.) sigmaZ=0.4e-3;
854 if(clRphi < 0.) clRphi=2.5e-3;
855 if(clZ < 0.) clZ=2.5e-5;
860 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
861 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
863 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
864 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
865 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
866 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
867 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
868 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
869 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
870 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
871 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
872 xyz[2]=sigmaRphi; // fSigmaY2
873 xyz[3]=sigmaZ; // fSigmaZ2
874 xyz[4]=tpcHit->fQ; // q
876 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
877 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
879 Int_t tracks[3]={tpcHit->Track(), -1, -1};
880 AliTPCcluster cluster(tracks,xyz);
882 clrow->InsertCluster(&cluster); nclusters++;
884 tpcHit = (AliTPChit*)NextHit();
887 } // end of loop over hits
889 } // end of loop over tracks
891 Int_t nrows=fTPCParam->GetNRow(isec);
892 for (Int_t irow=0; irow<nrows; irow++) {
893 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
894 if (!clrow) continue;
895 carray.StoreRow(isec,irow);
896 carray.ClearRow(isec,irow);
899 } // end of loop over sectors
901 cerr<<"Number of made clusters : "<<nclusters<<" \n";
903 carray.GetTree()->Write();
905 savedir->cd(); //switch back to the input file
909 //_________________________________________________________________
910 void AliTPC::Hits2ExactClustersSector(Int_t isec)
912 //--------------------------------------------------------
913 //calculate exact cross point of track and given pad row
914 //resulting values are expressed in "digit" coordinata
915 //--------------------------------------------------------
917 //-----------------------------------------------------------------
918 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
919 //-----------------------------------------------------------------
921 if (fClustersArray==0){
925 TParticle *particle; // pointer to a given particle
926 AliTPChit *tpcHit; // pointer to a sigle TPC hit
927 TClonesArray *particles; //pointer to the particle list
930 const Int_t kcmaxhits=30000;
931 TVector * xxxx = new TVector(kcmaxhits*4);
932 TVector & xxx = *xxxx;
933 Int_t maxhits = kcmaxhits;
934 //construct array for each padrow
935 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
936 fClustersArray->CreateRow(isec,i);
938 //---------------------------------------------------------------
939 // Get the access to the tracks
940 //---------------------------------------------------------------
942 TTree *tH = gAlice->TreeH();
943 Stat_t ntracks = tH->GetEntries();
944 particles=gAlice->Particles();
945 Int_t npart = particles->GetEntriesFast();
947 //------------------------------------------------------------
949 //------------------------------------------------------------
951 for(Int_t track=0;track<ntracks;track++){
955 // Get number of the TPC hits and a pointer
958 nhits=fHits->GetEntriesFast();
962 Int_t currentIndex=0;
963 Int_t lastrow=-1; //last writen row
964 for(Int_t hit=0;hit<nhits;hit++){
965 tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
966 if (tpcHit==0) continue;
967 sector=tpcHit->fSector; // sector number
968 if(sector != isec) continue;
969 ipart=tpcHit->Track();
970 if (ipart<npart) particle=(TParticle*)particles->UncheckedAt(ipart);
974 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
975 Int_t index[3]={1,isec,0};
976 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
977 if (currentrow<0) continue;
978 if (lastrow<0) lastrow=currentrow;
979 if (currentrow==lastrow){
980 if ( currentIndex>=maxhits){
982 xxx.ResizeTo(4*maxhits);
984 xxx(currentIndex*4)=x[0];
985 xxx(currentIndex*4+1)=x[1];
986 xxx(currentIndex*4+2)=x[2];
987 xxx(currentIndex*4+3)=tpcHit->fQ;
1003 for (Int_t index=0;index<currentIndex;index++){
1005 x=x2=x3=x4=xxx(index*4);
1013 sumy+=xxx(index*4+1);
1014 sumxy+=xxx(index*4+1)*x;
1015 sumx2y+=xxx(index*4+1)*x2;
1016 sumz+=xxx(index*4+2);
1017 sumxz+=xxx(index*4+2)*x;
1018 sumx2z+=xxx(index*4+2)*x2;
1019 sumq+=xxx(index*4+3);
1021 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1022 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1023 sumx2*(sumx*sumx3-sumx2*sumx2);
1025 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1026 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1027 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1028 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1030 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1031 sumx2*(sumx*sumx2y-sumx2*sumxy);
1032 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1033 sumx2*(sumx*sumx2z-sumx2*sumxz);
1035 Float_t y=detay/det+centralPad;
1036 Float_t z=detaz/det;
1037 Float_t by=detby/det; //y angle
1038 Float_t bz=detbz/det; //z angle
1039 sumy/=Float_t(currentIndex);
1040 sumz/=Float_t(currentIndex);
1041 AliComplexCluster cl;
1047 cl.fTracks[0]=ipart;
1049 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1050 if (row!=0) row->InsertCluster(&cl);
1053 } //end of calculating cluster for given row
1057 } // end of loop over hits
1058 } // end of loop over tracks
1059 //write padrows to tree
1060 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1061 fClustersArray->StoreRow(isec,ii);
1062 fClustersArray->ClearRow(isec,ii);
1068 //__________________________________________________________________
1069 void AliTPC::Hits2Digits()
1071 //----------------------------------------------------
1072 // Loop over all sectors
1073 //----------------------------------------------------
1076 printf("AliTPCParam MUST be created firstly\n");
1080 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1085 //_____________________________________________________________________________
1086 void AliTPC::Hits2DigitsSector(Int_t isec)
1088 //-------------------------------------------------------------------
1089 // TPC conversion from hits to digits.
1090 //-------------------------------------------------------------------
1092 //-----------------------------------------------------------------
1093 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1094 //-----------------------------------------------------------------
1096 //-------------------------------------------------------
1097 // Get the access to the track hits
1098 //-------------------------------------------------------
1101 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1102 Stat_t ntracks = tH->GetEntries();
1106 //-------------------------------------------
1107 // Only if there are any tracks...
1108 //-------------------------------------------
1112 printf("*** Processing sector number %d ***\n",isec);
1114 Int_t nrows =fTPCParam->GetNRow(isec);
1116 row= new TObjArray* [nrows];
1118 MakeSector(isec,nrows,tH,ntracks,row);
1120 //--------------------------------------------------------
1121 // Digitize this sector, row by row
1122 // row[i] is the pointer to the TObjArray of TVectors,
1123 // each one containing electrons accepted on this
1124 // row, assigned into tracks
1125 //--------------------------------------------------------
1129 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree();
1131 for (i=0;i<nrows;i++){
1133 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1135 DigitizeRow(i,isec,row);
1137 fDigitsArray->StoreRow(isec,i);
1139 Int_t ndig = dig->GetDigitSize();
1141 printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1143 fDigitsArray->ClearRow(isec,i);
1146 } // end of the sector digitization
1148 for(i=0;i<nrows;i++){
1153 delete [] row; // delete the array of pointers to TObjArray-s
1157 } // end of Hits2DigitsSector
1160 //_____________________________________________________________________________
1161 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1163 //-----------------------------------------------------------
1164 // Single row digitization, coupling from the neighbouring
1165 // rows taken into account
1166 //-----------------------------------------------------------
1168 //-----------------------------------------------------------------
1169 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1170 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1171 //-----------------------------------------------------------------
1174 Float_t zerosup = fTPCParam->GetZeroSup();
1175 Int_t nrows =fTPCParam->GetNRow(isec);
1176 fCurrentIndex[1]= isec;
1179 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1180 Int_t nofTbins = fTPCParam->GetMaxTBin();
1181 Int_t indexRange[4];
1183 // Integrated signal for this row
1184 // and a single track signal
1186 TMatrix *m1 = new TMatrix(0,nofPads,0,nofTbins); // integrated
1187 TMatrix *m2 = new TMatrix(0,nofPads,0,nofTbins); // single
1189 TMatrix &total = *m1;
1191 // Array of pointers to the label-signal list
1193 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1194 Float_t **pList = new Float_t* [nofDigits];
1198 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1202 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1203 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1204 for (Int_t row= row1;row<=row2;row++){
1205 Int_t nTracks= rows[row]->GetEntries();
1206 for (i1=0;i1<nTracks;i1++){
1207 fCurrentIndex[2]= row;
1208 fCurrentIndex[3]=irow;
1210 m2->Zero(); // clear single track signal matrix
1211 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1212 GetList(trackLabel,nofPads,m2,indexRange,pList);
1214 else GetSignal(rows[row],i1,0,m1,indexRange);
1220 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1221 for(Int_t ip=0;ip<nofPads;ip++){
1222 for(Int_t it=0;it<nofTbins;it++){
1224 Float_t q = total(ip,it);
1226 Int_t gi =it*nofPads+ip; // global index
1228 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
1232 if(q <=zerosup) continue; // do not fill zeros
1233 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1236 // "real" signal or electronic noise (list = -1)?
1239 for(Int_t j1=0;j1<3;j1++){
1240 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1245 <A NAME="AliDigits"></A>
1246 using of AliDigits object
1249 dig->SetDigitFast((Short_t)q,it,ip);
1250 if (fDigitsArray->IsSimulated())
1252 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1253 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1254 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1258 } // end of loop over time buckets
1259 } // end of lop over pads
1262 // This row has been digitized, delete nonused stuff
1265 for(lp=0;lp<nofDigits;lp++){
1266 if(pList[lp]) delete [] pList[lp];
1275 } // end of DigitizeRow
1277 //_____________________________________________________________________________
1279 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, TMatrix *m1, TMatrix *m2,
1283 //---------------------------------------------------------------
1284 // Calculates 2-D signal (pad,time) for a single track,
1285 // returns a pointer to the signal matrix and the track label
1286 // No digitization is performed at this level!!!
1287 //---------------------------------------------------------------
1289 //-----------------------------------------------------------------
1290 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1291 // Modified: Marian Ivanov
1292 //-----------------------------------------------------------------
1296 tv = (TVector*)p1->At(ntr); // pointer to a track
1299 Float_t label = v(0);
1300 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1302 Int_t nElectrons = (tv->GetNrows()-1)/4;
1303 indexRange[0]=9999; // min pad
1304 indexRange[1]=-1; // max pad
1305 indexRange[2]=9999; //min time
1306 indexRange[3]=-1; // max time
1308 // Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1310 TMatrix &signal = *m1;
1311 TMatrix &total = *m2;
1313 // Loop over all electrons
1315 for(Int_t nel=0; nel<nElectrons; nel++){
1317 Float_t aval = v(idx+4);
1318 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1319 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1320 Int_t n = fTPCParam->CalcResponse(xyz,fCurrentIndex,fCurrentIndex[3]);
1322 if (n>0) for (Int_t i =0; i<n; i++){
1323 Int_t *index = fTPCParam->GetResBin(i);
1324 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1325 if ( ( pad<(fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]))) && (pad>0)) {
1326 Int_t time=index[2];
1327 Float_t weight = fTPCParam->GetResWeight(i); //we normalise response to ADC channel
1328 weight *= eltoadcfac;
1330 if (m1!=0) signal(pad,time)+=weight;
1331 total(pad,time)+=weight;
1332 indexRange[0]=TMath::Min(indexRange[0],pad);
1333 indexRange[1]=TMath::Max(indexRange[1],pad);
1334 indexRange[2]=TMath::Min(indexRange[2],time);
1335 indexRange[3]=TMath::Max(indexRange[3],time);
1338 } // end of loop over electrons
1340 return label; // returns track label when finished
1343 //_____________________________________________________________________________
1344 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *indexRange,
1347 //----------------------------------------------------------------------
1348 // Updates the list of tracks contributing to digits for a given row
1349 //----------------------------------------------------------------------
1351 //-----------------------------------------------------------------
1352 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1353 //-----------------------------------------------------------------
1355 TMatrix &signal = *m;
1357 // lop over nonzero digits
1359 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1360 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1363 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1365 if(signal(ip,it)<0.5) continue;
1368 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1370 if(!pList[globalIndex]){
1373 // Create new list (6 elements - 3 signals and 3 labels),
1376 pList[globalIndex] = new Float_t [6];
1380 *pList[globalIndex] = -1.;
1381 *(pList[globalIndex]+1) = -1.;
1382 *(pList[globalIndex]+2) = -1.;
1383 *(pList[globalIndex]+3) = -1.;
1384 *(pList[globalIndex]+4) = -1.;
1385 *(pList[globalIndex]+5) = -1.;
1388 *pList[globalIndex] = label;
1389 *(pList[globalIndex]+3) = signal(ip,it);
1393 // check the signal magnitude
1395 Float_t highest = *(pList[globalIndex]+3);
1396 Float_t middle = *(pList[globalIndex]+4);
1397 Float_t lowest = *(pList[globalIndex]+5);
1400 // compare the new signal with already existing list
1403 if(signal(ip,it)<lowest) continue; // neglect this track
1407 if (signal(ip,it)>highest){
1408 *(pList[globalIndex]+5) = middle;
1409 *(pList[globalIndex]+4) = highest;
1410 *(pList[globalIndex]+3) = signal(ip,it);
1412 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1413 *(pList[globalIndex]+1) = *pList[globalIndex];
1414 *pList[globalIndex] = label;
1416 else if (signal(ip,it)>middle){
1417 *(pList[globalIndex]+5) = middle;
1418 *(pList[globalIndex]+4) = signal(ip,it);
1420 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1421 *(pList[globalIndex]+1) = label;
1424 *(pList[globalIndex]+5) = signal(ip,it);
1425 *(pList[globalIndex]+2) = label;
1429 } // end of loop over pads
1430 } // end of loop over time bins
1435 //___________________________________________________________________
1436 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1437 Stat_t ntracks,TObjArray **row)
1440 //-----------------------------------------------------------------
1441 // Prepares the sector digitization, creates the vectors of
1442 // tracks for each row of this sector. The track vector
1443 // contains the track label and the position of electrons.
1444 //-----------------------------------------------------------------
1446 //-----------------------------------------------------------------
1447 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1448 //-----------------------------------------------------------------
1450 Float_t gasgain = fTPCParam->GetGasGain();
1454 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1457 if (fHitType&2) branch = TH->GetBranch("TPC2");
1458 else branch = TH->GetBranch("TPC");
1461 //----------------------------------------------
1462 // Create TObjArray-s, one for each row,
1463 // each TObjArray will store the TVectors
1464 // of electrons, one TVector per each track.
1465 //----------------------------------------------
1467 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
1468 TVector **tracks = new TVector* [nrows]; //pointers to the track vectors
1469 for(i=0; i<nrows; i++){
1470 row[i] = new TObjArray;
1477 //--------------------------------------------------------------------
1478 // Loop over tracks, the "track" contains the full history
1479 //--------------------------------------------------------------------
1481 Int_t previousTrack,currentTrack;
1482 previousTrack = -1; // nothing to store so far!
1484 for(Int_t track=0;track<ntracks;track++){
1485 Bool_t isInSector=kTRUE;
1490 TBranch * br = TH->GetBranch("fTrackHitsInfo");
1491 br->GetEvent(track);
1492 AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
1493 for (UInt_t j=0;j<ar->GetSize();j++){
1494 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
1497 if (!isInSector) continue;
1499 branch->GetEntry(track); // get next track
1503 tpcHit = (AliTPChit*)FirstHit(-1);
1505 //--------------------------------------------------------------
1507 //--------------------------------------------------------------
1512 Int_t sector=tpcHit->fSector; // sector number
1513 // if(sector != isec) continue;
1515 tpcHit = (AliTPChit*) NextHit();
1519 currentTrack = tpcHit->Track(); // track number
1522 if(currentTrack != previousTrack){
1524 // store already filled fTrack
1526 for(i=0;i<nrows;i++){
1527 if(previousTrack != -1){
1528 if(nofElectrons[i]>0){
1529 TVector &v = *tracks[i];
1530 v(0) = previousTrack;
1531 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1532 row[i]->Add(tracks[i]);
1535 delete tracks[i]; // delete empty TVector
1541 tracks[i] = new TVector(481); // TVectors for the next fTrack
1543 } // end of loop over rows
1545 previousTrack=currentTrack; // update track label
1548 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1550 //---------------------------------------------------
1551 // Calculate the electron attachment probability
1552 //---------------------------------------------------
1555 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
1556 /fTPCParam->GetDriftV();
1558 Float_t attProb = fTPCParam->GetAttCoef()*
1559 fTPCParam->GetOxyCont()*time; // fraction!
1561 //-----------------------------------------------
1562 // Loop over electrons
1563 //-----------------------------------------------
1566 for(Int_t nel=0;nel<qI;nel++){
1567 // skip if electron lost due to the attachment
1568 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1572 xyz[3]= (Float_t) (-gasgain*TMath::Log(gRandom->Rndm()));
1575 TransportElectron(xyz,index); //MI change -august
1577 fTPCParam->GetPadRow(xyz,index); //MI change august
1578 rowNumber = index[2];
1579 //transform position to local digit coordinates
1580 //relative to nearest pad row
1581 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
1582 nofElectrons[rowNumber]++;
1583 //----------------------------------
1584 // Expand vector if necessary
1585 //----------------------------------
1586 if(nofElectrons[rowNumber]>120){
1587 Int_t range = tracks[rowNumber]->GetNrows();
1588 if((nofElectrons[rowNumber])>(range-1)/4){
1590 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
1594 TVector &v = *tracks[rowNumber];
1595 Int_t idx = 4*nofElectrons[rowNumber]-3;
1597 v(idx)= xyz[0]; // X - pad row coordinate
1598 v(idx+1)=xyz[1]; // Y - pad coordinate (along the pad-row)
1599 v(idx+2)=xyz[2]; // Z - time bin coordinate
1600 v(idx+3)=xyz[3]; // avalanche size
1601 } // end of loop over electrons
1603 tpcHit = (AliTPChit*)NextHit();
1605 } // end of loop over hits
1606 } // end of loop over tracks
1609 // store remaining track (the last one) if not empty
1612 for(i=0;i<nrows;i++){
1613 if(nofElectrons[i]>0){
1614 TVector &v = *tracks[i];
1615 v(0) = previousTrack;
1616 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1617 row[i]->Add(tracks[i]);
1626 delete [] nofElectrons;
1629 } // end of MakeSector
1632 //_____________________________________________________________________________
1636 // Initialise TPC detector after definition of geometry
1641 for(i=0;i<35;i++) printf("*");
1642 printf(" TPC_INIT ");
1643 for(i=0;i<35;i++) printf("*");
1646 for(i=0;i<80;i++) printf("*");
1650 //_____________________________________________________________________________
1651 void AliTPC::MakeBranch(Option_t* option)
1654 // Create Tree branches for the TPC.
1656 Int_t buffersize = 4000;
1657 char branchname[10];
1658 sprintf(branchname,"%s",GetName());
1660 AliDetector::MakeBranch(option);
1662 char *d = strstr(option,"D");
1664 if (fDigits && gAlice->TreeD() && d) {
1665 gAlice->TreeD()->Branch(branchname,&fDigits, buffersize);
1666 printf("Making Branch %s for digits\n",branchname);
1669 if (fHitType&2) MakeBranch2(option); // MI change 14.09.2000
1672 //_____________________________________________________________________________
1673 void AliTPC::ResetDigits()
1676 // Reset number of digits and the digits array for this detector
1679 if (fDigits) fDigits->Clear();
1682 //_____________________________________________________________________________
1683 void AliTPC::SetSecAL(Int_t sec)
1685 //---------------------------------------------------
1686 // Activate/deactivate selection for lower sectors
1687 //---------------------------------------------------
1689 //-----------------------------------------------------------------
1690 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1691 //-----------------------------------------------------------------
1696 //_____________________________________________________________________________
1697 void AliTPC::SetSecAU(Int_t sec)
1699 //----------------------------------------------------
1700 // Activate/deactivate selection for upper sectors
1701 //---------------------------------------------------
1703 //-----------------------------------------------------------------
1704 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1705 //-----------------------------------------------------------------
1710 //_____________________________________________________________________________
1711 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
1713 //----------------------------------------
1714 // Select active lower sectors
1715 //----------------------------------------
1717 //-----------------------------------------------------------------
1718 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1719 //-----------------------------------------------------------------
1729 //_____________________________________________________________________________
1730 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
1731 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
1732 Int_t s11 , Int_t s12)
1734 //--------------------------------
1735 // Select active upper sectors
1736 //--------------------------------
1738 //-----------------------------------------------------------------
1739 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1740 //-----------------------------------------------------------------
1756 //_____________________________________________________________________________
1757 void AliTPC::SetSens(Int_t sens)
1760 //-------------------------------------------------------------
1761 // Activates/deactivates the sensitive strips at the center of
1762 // the pad row -- this is for the space-point resolution calculations
1763 //-------------------------------------------------------------
1765 //-----------------------------------------------------------------
1766 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1767 //-----------------------------------------------------------------
1773 void AliTPC::SetSide(Float_t side=0.)
1775 // choice of the TPC side
1780 //____________________________________________________________________________
1781 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
1782 Float_t p2,Float_t p3)
1785 // gax mixture definition
1799 //_____________________________________________________________________________
1801 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
1804 // electron transport taking into account:
1806 // 2.ExB at the wires
1807 // 3. nonisochronity
1809 // xyz and index must be already transformed to system 1
1812 fTPCParam->Transform1to2(xyz,index);
1815 Float_t driftl=xyz[2];
1816 if(driftl<0.01) driftl=0.01;
1817 driftl=TMath::Sqrt(driftl);
1818 Float_t sigT = driftl*(fTPCParam->GetDiffT());
1819 Float_t sigL = driftl*(fTPCParam->GetDiffL());
1820 xyz[0]=gRandom->Gaus(xyz[0],sigT);
1821 xyz[1]=gRandom->Gaus(xyz[1],sigT);
1822 xyz[2]=gRandom->Gaus(xyz[2],sigL);
1826 if (fTPCParam->GetMWPCReadout()==kTRUE){
1828 fTPCParam->Transform2to2NearestWire(xyz,index);
1829 Float_t dx=xyz[0]-x1;
1830 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
1832 //add nonisochronity (not implemented yet)
1835 //_____________________________________________________________________________
1836 void AliTPC::Streamer(TBuffer &R__b)
1839 // Stream an object of class AliTPC.
1841 if (R__b.IsReading()) {
1842 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1843 AliDetector::Streamer(R__b);
1845 if (R__v < 2) return;
1850 R__b.WriteVersion(AliTPC::IsA());
1851 AliDetector::Streamer(R__b);
1858 ClassImp(AliTPCdigit)
1860 //_____________________________________________________________________________
1861 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
1865 // Creates a TPC digit object
1867 fSector = digits[0];
1868 fPadRow = digits[1];
1871 fSignal = digits[4];
1877 //_____________________________________________________________________________
1878 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
1882 // Creates a TPC hit object
1893 //________________________________________________________________________
1894 // Additional code because of the AliTPCTrackHits
1896 void AliTPC::MakeBranch2(Option_t *option)
1899 // Create a new branch in the current Root Tree
1900 // The branch of fHits is automatically split
1901 // MI change 14.09.2000
1902 if (fHitType&2==0) return;
1903 char branchname[10];
1904 sprintf(branchname,"%s2",GetName());
1906 // Get the pointer to the header
1907 char *cH = strstr(option,"H");
1909 if (fTrackHits && gAlice->TreeH() && cH) {
1910 // gAlice->TreeH()->Branch(branchname,&fTrackHits, fBufferSize);
1911 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHits,
1912 gAlice->TreeH(),fBufferSize,1);
1913 gAlice->TreeH()->GetListOfBranches()->Add(branch);
1914 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
1918 void AliTPC::SetTreeAddress()
1920 if (fHitType&1) AliDetector::SetTreeAddress();
1921 if (fHitType&2) SetTreeAddress2();
1924 void AliTPC::SetTreeAddress2()
1927 // Set branch address for the TrackHits Tree
1930 char branchname[20];
1931 sprintf(branchname,"%s2",GetName());
1933 // Branch address for hit tree
1934 TTree *treeH = gAlice->TreeH();
1936 branch = treeH->GetBranch(branchname);
1937 if (branch) branch->SetAddress(&fTrackHits);
1941 void AliTPC::FinishPrimary()
1943 if (fTrackHits) fTrackHits->FlushHitStack();
1947 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
1950 // add hit to the list
1951 TClonesArray &particles = *(gAlice->Particles());
1954 int primary = gAlice->GetPrimary(track);
1955 ((TParticle *)particles[primary])->SetBit(kKeepBit);
1959 gAlice->FlagTrack(track);
1961 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
1962 //if (hit->fTrack!=rtrack)
1963 // cout<<"bad track number\n";
1965 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
1966 hits[1],hits[2],(Int_t)hits[3]);
1969 void AliTPC::ResetHits()
1971 if (fHitType&1) AliDetector::ResetHits();
1972 if (fHitType&2) ResetHits2();
1975 void AliTPC::ResetHits2()
1979 if (fTrackHits) fTrackHits->Clear();
1982 AliHit* AliTPC::FirstHit(Int_t track)
1984 if (fHitType&2) return FirstHit2(track);
1985 return AliDetector::FirstHit(track);
1987 AliHit* AliTPC::NextHit()
1989 if (fHitType&2) return NextHit2();
1990 return AliDetector::NextHit();
1993 AliHit* AliTPC::FirstHit2(Int_t track)
1996 // Initialise the hit iterator
1997 // Return the address of the first hit for track
1998 // If track>=0 the track is read from disk
1999 // while if track<0 the first hit of the current
2000 // track is returned
2003 gAlice->ResetHits();
2004 gAlice->TreeH()->GetEvent(track);
2008 fTrackHits->First();
2009 return fTrackHits->GetHit();
2014 AliHit* AliTPC::NextHit2()
2017 //Return the next hit for the current track
2021 return fTrackHits->GetHit();
2027 void AliTPC::LoadPoints(Int_t)
2031 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2034 if(fHitType==1) AliDetector::LoadPoints(a);
2035 else LoadPoints2(a);
2042 void AliTPC::RemapTrackHitIDs(Int_t *map)
2044 if (!fTrackHits) return;
2045 AliObjectArray * arr = fTrackHits->fTrackHitsInfo;
2046 for (UInt_t i=0;i<arr->GetSize();i++){
2047 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2048 info->fTrackID = map[info->fTrackID];
2054 //_____________________________________________________________________________
2055 void AliTPC::LoadPoints2(Int_t)
2058 // Store x, y, z of all hits in memory
2060 if (fTrackHits == 0) return;
2062 Int_t nhits = fTrackHits->GetEntriesFast();
2063 if (nhits == 0) return;
2064 Int_t tracks = gAlice->GetNtrack();
2065 if (fPoints == 0) fPoints = new TObjArray(tracks);
2068 Int_t *ntrk=new Int_t[tracks];
2069 Int_t *limi=new Int_t[tracks];
2070 Float_t **coor=new Float_t*[tracks];
2071 for(Int_t i=0;i<tracks;i++) {
2077 AliPoints *points = 0;
2080 Int_t chunk=nhits/4+1;
2082 // Loop over all the hits and store their position
2084 ahit = FirstHit2(-1);
2085 //for (Int_t hit=0;hit<nhits;hit++) {
2087 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2088 trk=ahit->GetTrack();
2089 if(ntrk[trk]==limi[trk]) {
2091 // Initialise a new track
2092 fp=new Float_t[3*(limi[trk]+chunk)];
2094 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2095 delete [] coor[trk];
2102 fp[3*ntrk[trk] ] = ahit->X();
2103 fp[3*ntrk[trk]+1] = ahit->Y();
2104 fp[3*ntrk[trk]+2] = ahit->Z();
2109 for(trk=0; trk<tracks; ++trk) {
2111 points = new AliPoints();
2112 points->SetMarkerColor(GetMarkerColor());
2113 points->SetMarkerSize(GetMarkerSize());
2114 points->SetDetector(this);
2115 points->SetParticle(trk);
2116 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2117 fPoints->AddAt(points,trk);
2118 delete [] coor[trk];
2128 //_____________________________________________________________________________
2129 void AliTPC::LoadPoints3(Int_t)
2132 // Store x, y, z of all hits in memory
2133 // - only intersection point with pad row
2134 if (fTrackHits == 0) return;
2136 Int_t nhits = fTrackHits->GetEntriesFast();
2137 if (nhits == 0) return;
2138 Int_t tracks = gAlice->GetNtrack();
2139 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2140 fPoints->Expand(2*tracks);
2143 Int_t *ntrk=new Int_t[tracks];
2144 Int_t *limi=new Int_t[tracks];
2145 Float_t **coor=new Float_t*[tracks];
2146 for(Int_t i=0;i<tracks;i++) {
2152 AliPoints *points = 0;
2155 Int_t chunk=nhits/4+1;
2157 // Loop over all the hits and store their position
2159 ahit = FirstHit2(-1);
2160 //for (Int_t hit=0;hit<nhits;hit++) {
2164 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2165 trk=ahit->GetTrack();
2166 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2167 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2168 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2169 if (currentrow!=lastrow){
2170 lastrow = currentrow;
2171 //later calculate intersection point
2172 if(ntrk[trk]==limi[trk]) {
2174 // Initialise a new track
2175 fp=new Float_t[3*(limi[trk]+chunk)];
2177 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2178 delete [] coor[trk];
2185 fp[3*ntrk[trk] ] = ahit->X();
2186 fp[3*ntrk[trk]+1] = ahit->Y();
2187 fp[3*ntrk[trk]+2] = ahit->Z();
2194 for(trk=0; trk<tracks; ++trk) {
2196 points = new AliPoints();
2197 points->SetMarkerColor(GetMarkerColor()+1);
2198 points->SetMarkerStyle(5);
2199 points->SetMarkerSize(0.2);
2200 points->SetDetector(this);
2201 points->SetParticle(trk);
2202 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2203 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2204 fPoints->AddAt(points,tracks+trk);
2205 delete [] coor[trk];
2216 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2220 //fill clones array with intersection of current point with the
2225 const Int_t kcmaxhits=30000;
2226 TVector * xxxx = new TVector(kcmaxhits*4);
2227 TVector & xxx = *xxxx;
2228 Int_t maxhits = kcmaxhits;
2231 AliTPChit * tpcHit=0;
2232 tpcHit = (AliTPChit*)FirstHit2(-1);
2233 Int_t currentIndex=0;
2234 Int_t lastrow=-1; //last writen row
2237 if (tpcHit==0) continue;
2238 sector=tpcHit->fSector; // sector number
2239 ipart=tpcHit->Track();
2243 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2244 Int_t index[3]={1,sector,0};
2245 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2246 if (currentrow<0) continue;
2247 if (lastrow<0) lastrow=currentrow;
2248 if (currentrow==lastrow){
2249 if ( currentIndex>=maxhits){
2251 xxx.ResizeTo(4*maxhits);
2253 xxx(currentIndex*4)=x[0];
2254 xxx(currentIndex*4+1)=x[1];
2255 xxx(currentIndex*4+2)=x[2];
2256 xxx(currentIndex*4+3)=tpcHit->fQ;
2260 if (currentIndex>2){
2272 for (Int_t index=0;index<currentIndex;index++){
2274 x=x2=x3=x4=xxx(index*4);
2282 sumy+=xxx(index*4+1);
2283 sumxy+=xxx(index*4+1)*x;
2284 sumx2y+=xxx(index*4+1)*x2;
2285 sumz+=xxx(index*4+2);
2286 sumxz+=xxx(index*4+2)*x;
2287 sumx2z+=xxx(index*4+2)*x2;
2288 sumq+=xxx(index*4+3);
2290 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2291 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2292 sumx2*(sumx*sumx3-sumx2*sumx2);
2294 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2295 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2296 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2297 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2299 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2300 sumx2*(sumx*sumx2y-sumx2*sumxy);
2301 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2302 sumx2*(sumx*sumx2z-sumx2*sumxz);
2304 Float_t y=detay/det+centralPad;
2305 Float_t z=detaz/det;
2306 Float_t by=detby/det; //y angle
2307 Float_t bz=detbz/det; //z angle
2308 sumy/=Float_t(currentIndex);
2309 sumz/=Float_t(currentIndex);
2311 AliComplexCluster cl;
2317 cl.fTracks[0]=ipart;
2319 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2320 if (row!=0) row->InsertCluster(&cl);
2323 } //end of calculating cluster for given row
2325 } // end of loop over hits