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.30 2001/03/01 17:34:47 kowal2
19 Correction due to the accuracy problem
21 Revision 1.29 2001/02/28 16:34:40 kowal2
22 Protection against nonphysical values of the avalanche size,
25 Revision 1.28 2001/01/26 19:57:19 hristov
26 Major upgrade of AliRoot code
28 Revision 1.27 2001/01/13 17:29:33 kowal2
29 Sun compiler correction
31 Revision 1.26 2001/01/10 07:59:43 kowal2
32 Corrections to load points from the noncompressed hits.
34 Revision 1.25 2000/11/02 07:25:31 kowal2
35 Changes due to the new hit structure.
38 Revision 1.24 2000/10/05 16:06:09 kowal2
39 Forward declarations. Changes due to a new class AliComplexCluster.
41 Revision 1.23 2000/10/02 21:28:18 fca
42 Removal of useless dependecies via forward declarations
44 Revision 1.22 2000/07/10 20:57:39 hristov
45 Update of TPC code and macros by M.Kowalski
47 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
48 Changes to obey the coding rules
50 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
51 Splitted from AliTPCtracking
53 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
54 Changed parameter settings
56 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
58 Defaults loaded automatically (hard-wired)
59 Optional parameters can be set via macro called in the constructor
61 Revision 1.19 2000/04/18 19:00:59 fca
62 Small bug fixes to TPC files
64 Revision 1.18 2000/04/17 09:37:33 kowal2
65 removed obsolete AliTPCDigitsDisplay.C
67 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
69 New, experimental data structure from M. Ivanov
70 New tracking algorithm
71 Different pad geometry for different sectors
72 Digitization rewritten
74 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
75 Not used anymore - removed
77 Revision 1.17 2000/01/19 17:17:30 fca
78 Introducing a list of lists of hits -- more hits allowed for detector now
80 Revision 1.16 1999/11/05 09:29:23 fca
81 Accept only signals > 0
83 Revision 1.15 1999/10/08 06:26:53 fca
84 Removed ClustersIndex - not used anymore
86 Revision 1.14 1999/09/29 09:24:33 fca
87 Introduction of the Copyright and cvs Log
91 ///////////////////////////////////////////////////////////////////////////////
93 // Time Projection Chamber //
94 // This class contains the basic functions for the Time Projection Chamber //
95 // detector. Functions specific to one particular geometry are //
96 // contained in the derived classes //
100 <img src="picts/AliTPCClass.gif">
105 ///////////////////////////////////////////////////////////////////////////////
113 #include <TGeometry.h>
116 #include <TObjectTable.h>
117 #include "TParticle.h"
121 #include <iostream.h>
128 #include "AliTPCParamSR.h"
129 #include "AliTPCPRF2D.h"
130 #include "AliTPCRF1D.h"
131 #include "AliDigits.h"
132 #include "AliSimDigits.h"
133 #include "AliTPCTrackHits.h"
134 #include "AliPoints.h"
135 #include "AliArrayBranch.h"
138 #include "AliTPCDigitsArray.h"
139 #include "AliComplexCluster.h"
140 #include "AliClusters.h"
141 #include "AliTPCClustersRow.h"
142 #include "AliTPCClustersArray.h"
144 #include "AliTPCcluster.h"
145 #include "AliTPCclusterer.h"
146 #include "AliTPCtracker.h"
148 #include <TInterpreter.h>
155 //_____________________________________________________________________________
159 // Default constructor
174 //_____________________________________________________________________________
175 AliTPC::AliTPC(const char *name, const char *title)
176 : AliDetector(name,title)
179 // Standard constructor
183 // Initialise arrays of hits and digits
184 fHits = new TClonesArray("AliTPChit", 176);
185 gAlice->AddHitList(fHits);
190 fTrackHits = new AliTPCTrackHits; //MI - 13.09.2000
191 fTrackHits->SetHitPrecision(0.002);
192 fTrackHits->SetStepPrecision(0.003);
193 fTrackHits->SetMaxDistance(100);
197 // Initialise counters
203 // Initialise color attributes
204 SetMarkerColor(kYellow);
207 // Set TPC parameters
210 if (!strcmp(title,"Default")) {
211 fTPCParam = new AliTPCParamSR;
213 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
219 //_____________________________________________________________________________
230 delete fTrackHits; //MI 15.09.2000
233 //_____________________________________________________________________________
234 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
237 // Add a hit to the list
239 // TClonesArray &lhits = *fHits;
240 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
242 TClonesArray &lhits = *fHits;
243 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
246 AddHit2(track,vol,hits);
249 //_____________________________________________________________________________
250 void AliTPC::BuildGeometry()
254 // Build TPC ROOT TNode geometry for the event display
259 const int kColorTPC=19;
260 char name[5], title[25];
261 const Double_t kDegrad=TMath::Pi()/180;
262 const Double_t kRaddeg=180./TMath::Pi();
265 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
266 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
268 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
269 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
271 Int_t nLo = fTPCParam->GetNInnerSector()/2;
272 Int_t nHi = fTPCParam->GetNOuterSector()/2;
274 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
275 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
276 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
277 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
280 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
281 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
287 // Get ALICE top node
290 nTop=gAlice->GetGeometry()->GetNode("alice");
294 rl = fTPCParam->GetInnerRadiusLow();
295 ru = fTPCParam->GetInnerRadiusUp();
299 sprintf(name,"LS%2.2d",i);
301 sprintf(title,"TPC low sector %3d",i);
304 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
305 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
306 tubs->SetNumberOfDivisions(1);
308 nNode = new TNode(name,title,name,0,0,0,"");
309 nNode->SetLineColor(kColorTPC);
315 rl = fTPCParam->GetOuterRadiusLow();
316 ru = fTPCParam->GetOuterRadiusUp();
319 sprintf(name,"US%2.2d",i);
321 sprintf(title,"TPC upper sector %d",i);
323 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
324 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
325 tubs->SetNumberOfDivisions(1);
327 nNode = new TNode(name,title,name,0,0,0,"");
328 nNode->SetLineColor(kColorTPC);
334 //_____________________________________________________________________________
335 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
338 // Calculate distance from TPC to mouse on the display
344 void AliTPC::Clusters2Tracks(TFile *of) {
345 //-----------------------------------------------------------------
346 // This is a track finder.
347 //-----------------------------------------------------------------
348 AliTPCtracker::Clusters2Tracks(fTPCParam,of);
351 //_____________________________________________________________________________
352 void AliTPC::CreateMaterials()
354 //-----------------------------------------------
355 // Create Materials for for TPC simulations
356 //-----------------------------------------------
358 //-----------------------------------------------------------------
359 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
360 //-----------------------------------------------------------------
362 Int_t iSXFLD=gAlice->Field()->Integ();
363 Float_t sXMGMX=gAlice->Field()->Max();
365 Float_t amat[5]; // atomic numbers
366 Float_t zmat[5]; // z
367 Float_t wmat[5]; // proportions
373 //***************** Gases *************************
375 //-------------------------------------------------
377 //-------------------------------------------------
388 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
398 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
401 //--------------------------------------------------------------
403 //--------------------------------------------------------------
420 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
422 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
437 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
439 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
455 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
457 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
459 //----------------------------------------------------------------
460 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
461 //----------------------------------------------------------------
467 Float_t rho,absl,X0,buf[1];
471 for(nc = 0;nc<fNoComp;nc++)
474 // retrive material constants
476 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
481 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
483 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
484 density += fMixtProp[nc]*rho; // density of the mixture
488 // mixture proportions by weight!
490 for(nc = 0;nc<fNoComp;nc++)
493 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
495 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
496 apure[nnc] : amol[nnc])/am;
500 // Drift gases 1 - nonsensitive, 2 - sensitive
502 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
503 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
512 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
515 //----------------------------------------------------------------------
517 //----------------------------------------------------------------------
539 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
561 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
579 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
597 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
599 // G10 60% SiO2 + 40% epoxy, I use A and Z for SiO2
612 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz
614 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
616 AliMaterial(39,"G10",amat[0],zmat[0],density,999.,999.);
625 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
634 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
643 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
661 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
680 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
684 //----------------------------------------------------------
685 // tracking media for gases
686 //----------------------------------------------------------
688 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
689 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
690 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
691 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
693 //-----------------------------------------------------------
694 // tracking media for solids
695 //-----------------------------------------------------------
697 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
698 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
699 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
700 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
701 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
702 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
703 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
704 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
705 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
706 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
711 void AliTPC::Digits2Clusters(TFile *of)
713 //-----------------------------------------------------------------
714 // This is a simple cluster finder.
715 //-----------------------------------------------------------------
716 AliTPCclusterer::Digits2Clusters(fTPCParam,of);
719 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
720 extern Double_t SigmaZ2(Double_t, Double_t);
721 //_____________________________________________________________________________
722 void AliTPC::Hits2Clusters(TFile *of)
724 //--------------------------------------------------------
725 // TPC simple cluster generator from hits
726 // obtained from the TPC Fast Simulator
727 // The point errors are taken from the parametrization
728 //--------------------------------------------------------
730 //-----------------------------------------------------------------
731 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
732 //-----------------------------------------------------------------
733 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
734 // Jouri.Belikov@cern.ch
735 //----------------------------------------------------------------
737 /////////////////////////////////////////////////////////////////////////////
739 //---------------------------------------------------------------------
740 // ALICE TPC Cluster Parameters
741 //--------------------------------------------------------------------
745 // Cluster width in rphi
746 const Float_t kACrphi=0.18322;
747 const Float_t kBCrphi=0.59551e-3;
748 const Float_t kCCrphi=0.60952e-1;
749 // Cluster width in z
750 const Float_t kACz=0.19081;
751 const Float_t kBCz=0.55938e-3;
752 const Float_t kCCz=0.30428;
754 TDirectory *savedir=gDirectory;
757 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
762 printf("AliTPCParam MUST be created firstly\n");
766 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
768 TParticle *particle; // pointer to a given particle
769 AliTPChit *tpcHit; // pointer to a sigle TPC hit
773 Float_t pl,pt,tanth,rpad,ratio;
776 //---------------------------------------------------------------
777 // Get the access to the tracks
778 //---------------------------------------------------------------
780 TTree *tH = gAlice->TreeH();
781 Stat_t ntracks = tH->GetEntries();
783 //Switch to the output file
786 fTPCParam->Write(fTPCParam->GetTitle());
787 AliTPCClustersArray carray;
788 carray.Setup(fTPCParam);
789 carray.SetClusterType("AliTPCcluster");
792 Int_t nclusters=0; //cluster counter
794 //------------------------------------------------------------
795 // Loop over all sectors (72 sectors for 20 deg
796 // segmentation for both lower and upper sectors)
797 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
798 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
800 // First cluster for sector 0 starts at "0"
801 //------------------------------------------------------------
803 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
805 fTPCParam->AdjustCosSin(isec,cph,sph);
807 //------------------------------------------------------------
809 //------------------------------------------------------------
811 for(Int_t track=0;track<ntracks;track++){
815 // Get number of the TPC hits
817 // nhits=fHits->GetEntriesFast();
820 tpcHit = (AliTPChit*)FirstHit(-1);
824 // for(Int_t hit=0;hit<nhits;hit++){
825 //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
829 if (tpcHit->fQ == 0.) {
830 tpcHit = (AliTPChit*) NextHit();
831 continue; //information about track (I.Belikov)
833 sector=tpcHit->fSector; // sector number
836 // if(sector != isec) continue; //terminate iteration
839 tpcHit = (AliTPChit*) NextHit();
842 ipart=tpcHit->Track();
843 particle=gAlice->Particle(ipart);
846 if(pt < 1.e-9) pt=1.e-9;
848 tanth = TMath::Abs(tanth);
849 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
850 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
852 // space-point resolutions
854 sigmaRphi=SigmaY2(rpad,tanth,pt);
855 sigmaZ =SigmaZ2(rpad,tanth );
859 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
860 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
862 // temporary protection
864 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
865 if(sigmaZ < 0.) sigmaZ=0.4e-3;
866 if(clRphi < 0.) clRphi=2.5e-3;
867 if(clZ < 0.) clZ=2.5e-5;
872 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
873 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
875 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
876 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
877 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
878 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
879 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
880 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
881 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
882 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
883 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
884 xyz[2]=sigmaRphi; // fSigmaY2
885 xyz[3]=sigmaZ; // fSigmaZ2
886 xyz[4]=tpcHit->fQ; // q
888 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
889 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
891 Int_t tracks[3]={tpcHit->Track(), -1, -1};
892 AliTPCcluster cluster(tracks,xyz);
894 clrow->InsertCluster(&cluster); nclusters++;
896 tpcHit = (AliTPChit*)NextHit();
899 } // end of loop over hits
901 } // end of loop over tracks
903 Int_t nrows=fTPCParam->GetNRow(isec);
904 for (Int_t irow=0; irow<nrows; irow++) {
905 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
906 if (!clrow) continue;
907 carray.StoreRow(isec,irow);
908 carray.ClearRow(isec,irow);
911 } // end of loop over sectors
913 cerr<<"Number of made clusters : "<<nclusters<<" \n";
915 carray.GetTree()->Write();
917 savedir->cd(); //switch back to the input file
921 //_________________________________________________________________
922 void AliTPC::Hits2ExactClustersSector(Int_t isec)
924 //--------------------------------------------------------
925 //calculate exact cross point of track and given pad row
926 //resulting values are expressed in "digit" coordinata
927 //--------------------------------------------------------
929 //-----------------------------------------------------------------
930 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
931 //-----------------------------------------------------------------
933 if (fClustersArray==0){
937 TParticle *particle; // pointer to a given particle
938 AliTPChit *tpcHit; // pointer to a sigle TPC hit
941 const Int_t kcmaxhits=30000;
942 TVector * xxxx = new TVector(kcmaxhits*4);
943 TVector & xxx = *xxxx;
944 Int_t maxhits = kcmaxhits;
945 //construct array for each padrow
946 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
947 fClustersArray->CreateRow(isec,i);
949 //---------------------------------------------------------------
950 // Get the access to the tracks
951 //---------------------------------------------------------------
953 TTree *tH = gAlice->TreeH();
954 Stat_t ntracks = tH->GetEntries();
955 Int_t npart = gAlice->GetNtrack();
957 //------------------------------------------------------------
959 //------------------------------------------------------------
961 for(Int_t track=0;track<ntracks;track++){
965 // Get number of the TPC hits and a pointer
968 nhits=fHits->GetEntriesFast();
972 Int_t currentIndex=0;
973 Int_t lastrow=-1; //last writen row
974 for(Int_t hit=0;hit<nhits;hit++){
975 tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
976 if (tpcHit==0) continue;
977 sector=tpcHit->fSector; // sector number
978 if(sector != isec) continue;
979 ipart=tpcHit->Track();
980 if (ipart<npart) particle=gAlice->Particle(ipart);
984 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
985 Int_t index[3]={1,isec,0};
986 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
987 if (currentrow<0) continue;
988 if (lastrow<0) lastrow=currentrow;
989 if (currentrow==lastrow){
990 if ( currentIndex>=maxhits){
992 xxx.ResizeTo(4*maxhits);
994 xxx(currentIndex*4)=x[0];
995 xxx(currentIndex*4+1)=x[1];
996 xxx(currentIndex*4+2)=x[2];
997 xxx(currentIndex*4+3)=tpcHit->fQ;
1001 if (currentIndex>2){
1013 for (Int_t index=0;index<currentIndex;index++){
1015 x=x2=x3=x4=xxx(index*4);
1023 sumy+=xxx(index*4+1);
1024 sumxy+=xxx(index*4+1)*x;
1025 sumx2y+=xxx(index*4+1)*x2;
1026 sumz+=xxx(index*4+2);
1027 sumxz+=xxx(index*4+2)*x;
1028 sumx2z+=xxx(index*4+2)*x2;
1029 sumq+=xxx(index*4+3);
1031 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1032 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1033 sumx2*(sumx*sumx3-sumx2*sumx2);
1035 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1036 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1037 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1038 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1040 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1041 sumx2*(sumx*sumx2y-sumx2*sumxy);
1042 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1043 sumx2*(sumx*sumx2z-sumx2*sumxz);
1045 Float_t y=detay/det+centralPad;
1046 Float_t z=detaz/det;
1047 Float_t by=detby/det; //y angle
1048 Float_t bz=detbz/det; //z angle
1049 sumy/=Float_t(currentIndex);
1050 sumz/=Float_t(currentIndex);
1051 AliComplexCluster cl;
1057 cl.fTracks[0]=ipart;
1059 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1060 if (row!=0) row->InsertCluster(&cl);
1063 } //end of calculating cluster for given row
1067 } // end of loop over hits
1068 } // end of loop over tracks
1069 //write padrows to tree
1070 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1071 fClustersArray->StoreRow(isec,ii);
1072 fClustersArray->ClearRow(isec,ii);
1077 //___________________________________________
1078 void AliTPC::SDigits2Digits()
1080 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1081 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1082 AliTPCPRF2D * prfouter = new AliTPCPRF2D;
1083 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1085 TDirectory *cwd = gDirectory;
1086 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1087 rf->SetOffset(3*param->GetZSigma());
1089 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1091 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n";
1094 prfinner->Read("prf_07504_Gati_056068_d02");
1095 prfouter->Read("prf_10006_Gati_047051_d03");
1099 param->SetInnerPRF(prfinner);
1100 param->SetOuterPRF(prfouter);
1101 param->SetTimeRF(rf);
1105 cerr<<"Digitizing TPC...\n";
1107 //setup TPCDigitsArray
1108 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1109 arr->SetClass("AliSimDigits");
1113 arr->MakeTree(fDigitsFile);
1115 SetDigitsArray(arr);
1119 // Hits2DigitsSector(1);
1120 // Hits2DigitsSector(2);
1121 // Hits2DigitsSector(3);
1122 // Hits2DigitsSector(1+18);
1123 // Hits2DigitsSector(2+18);
1124 // Hits2DigitsSector(3+18);
1126 // Hits2DigitsSector(36+1);
1127 // Hits2DigitsSector(36+2);
1128 // Hits2DigitsSector(36+3);
1129 // Hits2DigitsSector(36+1+18);
1130 // Hits2DigitsSector(36+2+18);
1131 // Hits2DigitsSector(36+3+18);
1136 sprintf(treeName,"TreeD_%s",param->GetTitle());
1137 GetDigitsArray()->GetTree()->Write(treeName,TObject::kOverwrite);
1140 //__________________________________________________________________
1141 void AliTPC::Hits2Digits()
1143 //----------------------------------------------------
1144 // Loop over all sectors
1145 //----------------------------------------------------
1148 printf("AliTPCParam MUST be created firstly\n");
1152 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1157 //_____________________________________________________________________________
1158 void AliTPC::Hits2DigitsSector(Int_t isec)
1160 //-------------------------------------------------------------------
1161 // TPC conversion from hits to digits.
1162 //-------------------------------------------------------------------
1164 //-----------------------------------------------------------------
1165 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1166 //-----------------------------------------------------------------
1168 //-------------------------------------------------------
1169 // Get the access to the track hits
1170 //-------------------------------------------------------
1173 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1174 Stat_t ntracks = tH->GetEntries();
1178 //-------------------------------------------
1179 // Only if there are any tracks...
1180 //-------------------------------------------
1184 printf("*** Processing sector number %d ***\n",isec);
1186 Int_t nrows =fTPCParam->GetNRow(isec);
1188 row= new TObjArray* [nrows];
1190 MakeSector(isec,nrows,tH,ntracks,row);
1192 //--------------------------------------------------------
1193 // Digitize this sector, row by row
1194 // row[i] is the pointer to the TObjArray of TVectors,
1195 // each one containing electrons accepted on this
1196 // row, assigned into tracks
1197 //--------------------------------------------------------
1201 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1203 for (i=0;i<nrows;i++){
1205 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1207 DigitizeRow(i,isec,row);
1209 fDigitsArray->StoreRow(isec,i);
1211 Int_t ndig = dig->GetDigitSize();
1213 printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1215 fDigitsArray->ClearRow(isec,i);
1218 } // end of the sector digitization
1220 for(i=0;i<nrows;i++){
1225 delete [] row; // delete the array of pointers to TObjArray-s
1229 } // end of Hits2DigitsSector
1232 //_____________________________________________________________________________
1233 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1235 //-----------------------------------------------------------
1236 // Single row digitization, coupling from the neighbouring
1237 // rows taken into account
1238 //-----------------------------------------------------------
1240 //-----------------------------------------------------------------
1241 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1242 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1243 //-----------------------------------------------------------------
1246 Float_t zerosup = fTPCParam->GetZeroSup();
1247 Int_t nrows =fTPCParam->GetNRow(isec);
1248 fCurrentIndex[1]= isec;
1251 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1252 Int_t nofTbins = fTPCParam->GetMaxTBin();
1253 Int_t indexRange[4];
1255 // Integrated signal for this row
1256 // and a single track signal
1258 TMatrix *m1 = new TMatrix(0,nofPads,0,nofTbins); // integrated
1259 TMatrix *m2 = new TMatrix(0,nofPads,0,nofTbins); // single
1261 TMatrix &total = *m1;
1263 // Array of pointers to the label-signal list
1265 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1266 Float_t **pList = new Float_t* [nofDigits];
1270 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1274 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1275 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1276 for (Int_t row= row1;row<=row2;row++){
1277 Int_t nTracks= rows[row]->GetEntries();
1278 for (i1=0;i1<nTracks;i1++){
1279 fCurrentIndex[2]= row;
1280 fCurrentIndex[3]=irow;
1282 m2->Zero(); // clear single track signal matrix
1283 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1284 GetList(trackLabel,nofPads,m2,indexRange,pList);
1286 else GetSignal(rows[row],i1,0,m1,indexRange);
1292 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1293 for(Int_t ip=0;ip<nofPads;ip++){
1294 for(Int_t it=0;it<nofTbins;it++){
1296 Float_t q = total(ip,it);
1298 Int_t gi =it*nofPads+ip; // global index
1300 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
1304 if(q <=zerosup) continue; // do not fill zeros
1305 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1308 // "real" signal or electronic noise (list = -1)?
1311 for(Int_t j1=0;j1<3;j1++){
1312 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1317 <A NAME="AliDigits"></A>
1318 using of AliDigits object
1321 dig->SetDigitFast((Short_t)q,it,ip);
1322 if (fDigitsArray->IsSimulated())
1324 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1325 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1326 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1330 } // end of loop over time buckets
1331 } // end of lop over pads
1334 // This row has been digitized, delete nonused stuff
1337 for(lp=0;lp<nofDigits;lp++){
1338 if(pList[lp]) delete [] pList[lp];
1347 } // end of DigitizeRow
1349 //_____________________________________________________________________________
1351 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, TMatrix *m1, TMatrix *m2,
1355 //---------------------------------------------------------------
1356 // Calculates 2-D signal (pad,time) for a single track,
1357 // returns a pointer to the signal matrix and the track label
1358 // No digitization is performed at this level!!!
1359 //---------------------------------------------------------------
1361 //-----------------------------------------------------------------
1362 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1363 // Modified: Marian Ivanov
1364 //-----------------------------------------------------------------
1368 tv = (TVector*)p1->At(ntr); // pointer to a track
1371 Float_t label = v(0);
1372 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1374 Int_t nElectrons = (tv->GetNrows()-1)/4;
1375 indexRange[0]=9999; // min pad
1376 indexRange[1]=-1; // max pad
1377 indexRange[2]=9999; //min time
1378 indexRange[3]=-1; // max time
1380 // Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1382 TMatrix &signal = *m1;
1383 TMatrix &total = *m2;
1385 // Loop over all electrons
1387 for(Int_t nel=0; nel<nElectrons; nel++){
1389 Float_t aval = v(idx+4);
1390 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1391 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1392 Int_t n = fTPCParam->CalcResponse(xyz,fCurrentIndex,fCurrentIndex[3]);
1394 if (n>0) for (Int_t i =0; i<n; i++){
1395 Int_t *index = fTPCParam->GetResBin(i);
1396 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1397 if ( ( pad<(fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]))) && (pad>0)) {
1398 Int_t time=index[2];
1399 Float_t weight = fTPCParam->GetResWeight(i); //we normalise response to ADC channel
1400 weight *= eltoadcfac;
1402 if (m1!=0) signal(pad,time)+=weight;
1403 total(pad,time)+=weight;
1404 indexRange[0]=TMath::Min(indexRange[0],pad);
1405 indexRange[1]=TMath::Max(indexRange[1],pad);
1406 indexRange[2]=TMath::Min(indexRange[2],time);
1407 indexRange[3]=TMath::Max(indexRange[3],time);
1410 } // end of loop over electrons
1412 return label; // returns track label when finished
1415 //_____________________________________________________________________________
1416 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *indexRange,
1419 //----------------------------------------------------------------------
1420 // Updates the list of tracks contributing to digits for a given row
1421 //----------------------------------------------------------------------
1423 //-----------------------------------------------------------------
1424 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1425 //-----------------------------------------------------------------
1427 TMatrix &signal = *m;
1429 // lop over nonzero digits
1431 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1432 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1435 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1437 if(signal(ip,it)<0.5) continue;
1440 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1442 if(!pList[globalIndex]){
1445 // Create new list (6 elements - 3 signals and 3 labels),
1448 pList[globalIndex] = new Float_t [6];
1452 *pList[globalIndex] = -1.;
1453 *(pList[globalIndex]+1) = -1.;
1454 *(pList[globalIndex]+2) = -1.;
1455 *(pList[globalIndex]+3) = -1.;
1456 *(pList[globalIndex]+4) = -1.;
1457 *(pList[globalIndex]+5) = -1.;
1460 *pList[globalIndex] = label;
1461 *(pList[globalIndex]+3) = signal(ip,it);
1465 // check the signal magnitude
1467 Float_t highest = *(pList[globalIndex]+3);
1468 Float_t middle = *(pList[globalIndex]+4);
1469 Float_t lowest = *(pList[globalIndex]+5);
1472 // compare the new signal with already existing list
1475 if(signal(ip,it)<lowest) continue; // neglect this track
1479 if (signal(ip,it)>highest){
1480 *(pList[globalIndex]+5) = middle;
1481 *(pList[globalIndex]+4) = highest;
1482 *(pList[globalIndex]+3) = signal(ip,it);
1484 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1485 *(pList[globalIndex]+1) = *pList[globalIndex];
1486 *pList[globalIndex] = label;
1488 else if (signal(ip,it)>middle){
1489 *(pList[globalIndex]+5) = middle;
1490 *(pList[globalIndex]+4) = signal(ip,it);
1492 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1493 *(pList[globalIndex]+1) = label;
1496 *(pList[globalIndex]+5) = signal(ip,it);
1497 *(pList[globalIndex]+2) = label;
1501 } // end of loop over pads
1502 } // end of loop over time bins
1507 //___________________________________________________________________
1508 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1509 Stat_t ntracks,TObjArray **row)
1512 //-----------------------------------------------------------------
1513 // Prepares the sector digitization, creates the vectors of
1514 // tracks for each row of this sector. The track vector
1515 // contains the track label and the position of electrons.
1516 //-----------------------------------------------------------------
1518 //-----------------------------------------------------------------
1519 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1520 //-----------------------------------------------------------------
1522 Float_t gasgain = fTPCParam->GetGasGain();
1526 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1529 if (fHitType&2) branch = TH->GetBranch("TPC2");
1530 else branch = TH->GetBranch("TPC");
1533 //----------------------------------------------
1534 // Create TObjArray-s, one for each row,
1535 // each TObjArray will store the TVectors
1536 // of electrons, one TVector per each track.
1537 //----------------------------------------------
1539 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
1540 TVector **tracks = new TVector* [nrows]; //pointers to the track vectors
1541 for(i=0; i<nrows; i++){
1542 row[i] = new TObjArray;
1549 //--------------------------------------------------------------------
1550 // Loop over tracks, the "track" contains the full history
1551 //--------------------------------------------------------------------
1553 Int_t previousTrack,currentTrack;
1554 previousTrack = -1; // nothing to store so far!
1556 for(Int_t track=0;track<ntracks;track++){
1557 Bool_t isInSector=kTRUE;
1562 TBranch * br = TH->GetBranch("fTrackHitsInfo");
1563 br->GetEvent(track);
1564 AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
1565 for (UInt_t j=0;j<ar->GetSize();j++){
1566 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
1569 if (!isInSector) continue;
1571 branch->GetEntry(track); // get next track
1575 tpcHit = (AliTPChit*)FirstHit(-1);
1577 //--------------------------------------------------------------
1579 //--------------------------------------------------------------
1584 Int_t sector=tpcHit->fSector; // sector number
1585 // if(sector != isec) continue;
1587 tpcHit = (AliTPChit*) NextHit();
1591 currentTrack = tpcHit->Track(); // track number
1594 if(currentTrack != previousTrack){
1596 // store already filled fTrack
1598 for(i=0;i<nrows;i++){
1599 if(previousTrack != -1){
1600 if(nofElectrons[i]>0){
1601 TVector &v = *tracks[i];
1602 v(0) = previousTrack;
1603 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1604 row[i]->Add(tracks[i]);
1607 delete tracks[i]; // delete empty TVector
1613 tracks[i] = new TVector(481); // TVectors for the next fTrack
1615 } // end of loop over rows
1617 previousTrack=currentTrack; // update track label
1620 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1622 //---------------------------------------------------
1623 // Calculate the electron attachment probability
1624 //---------------------------------------------------
1627 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
1628 /fTPCParam->GetDriftV();
1630 Float_t attProb = fTPCParam->GetAttCoef()*
1631 fTPCParam->GetOxyCont()*time; // fraction!
1633 //-----------------------------------------------
1634 // Loop over electrons
1635 //-----------------------------------------------
1638 for(Int_t nel=0;nel<qI;nel++){
1639 // skip if electron lost due to the attachment
1640 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1645 // protection for the nonphysical avalanche size (10**6 maximum)
1647 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
1648 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
1651 TransportElectron(xyz,index); //MI change -august
1653 fTPCParam->GetPadRow(xyz,index); //MI change august
1654 rowNumber = index[2];
1655 //transform position to local digit coordinates
1656 //relative to nearest pad row
1657 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
1658 nofElectrons[rowNumber]++;
1659 //----------------------------------
1660 // Expand vector if necessary
1661 //----------------------------------
1662 if(nofElectrons[rowNumber]>120){
1663 Int_t range = tracks[rowNumber]->GetNrows();
1664 if((nofElectrons[rowNumber])>(range-1)/4){
1666 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
1670 TVector &v = *tracks[rowNumber];
1671 Int_t idx = 4*nofElectrons[rowNumber]-3;
1673 v(idx)= xyz[0]; // X - pad row coordinate
1674 v(idx+1)=xyz[1]; // Y - pad coordinate (along the pad-row)
1675 v(idx+2)=xyz[2]; // Z - time bin coordinate
1676 v(idx+3)=xyz[3]; // avalanche size
1677 } // end of loop over electrons
1679 tpcHit = (AliTPChit*)NextHit();
1681 } // end of loop over hits
1682 } // end of loop over tracks
1685 // store remaining track (the last one) if not empty
1688 for(i=0;i<nrows;i++){
1689 if(nofElectrons[i]>0){
1690 TVector &v = *tracks[i];
1691 v(0) = previousTrack;
1692 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1693 row[i]->Add(tracks[i]);
1702 delete [] nofElectrons;
1705 } // end of MakeSector
1708 //_____________________________________________________________________________
1712 // Initialise TPC detector after definition of geometry
1717 for(i=0;i<35;i++) printf("*");
1718 printf(" TPC_INIT ");
1719 for(i=0;i<35;i++) printf("*");
1722 for(i=0;i<80;i++) printf("*");
1726 //_____________________________________________________________________________
1727 void AliTPC::MakeBranch(Option_t* option, char *file)
1730 // Create Tree branches for the TPC.
1732 Int_t buffersize = 4000;
1733 char branchname[10];
1734 sprintf(branchname,"%s",GetName());
1736 AliDetector::MakeBranch(option,file);
1738 char *d = strstr(option,"D");
1740 if (fDigits && gAlice->TreeD() && d) {
1741 gAlice->MakeBranchInTree(gAlice->TreeD(),
1742 branchname, &fDigits, buffersize, file) ;
1745 if (fHitType&2) MakeBranch2(option,file); // MI change 14.09.2000
1748 //_____________________________________________________________________________
1749 void AliTPC::ResetDigits()
1752 // Reset number of digits and the digits array for this detector
1755 if (fDigits) fDigits->Clear();
1758 //_____________________________________________________________________________
1759 void AliTPC::SetSecAL(Int_t sec)
1761 //---------------------------------------------------
1762 // Activate/deactivate selection for lower sectors
1763 //---------------------------------------------------
1765 //-----------------------------------------------------------------
1766 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1767 //-----------------------------------------------------------------
1772 //_____________________________________________________________________________
1773 void AliTPC::SetSecAU(Int_t sec)
1775 //----------------------------------------------------
1776 // Activate/deactivate selection for upper sectors
1777 //---------------------------------------------------
1779 //-----------------------------------------------------------------
1780 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1781 //-----------------------------------------------------------------
1786 //_____________________________________________________________________________
1787 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
1789 //----------------------------------------
1790 // Select active lower sectors
1791 //----------------------------------------
1793 //-----------------------------------------------------------------
1794 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1795 //-----------------------------------------------------------------
1805 //_____________________________________________________________________________
1806 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
1807 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
1808 Int_t s11 , Int_t s12)
1810 //--------------------------------
1811 // Select active upper sectors
1812 //--------------------------------
1814 //-----------------------------------------------------------------
1815 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1816 //-----------------------------------------------------------------
1832 //_____________________________________________________________________________
1833 void AliTPC::SetSens(Int_t sens)
1836 //-------------------------------------------------------------
1837 // Activates/deactivates the sensitive strips at the center of
1838 // the pad row -- this is for the space-point resolution calculations
1839 //-------------------------------------------------------------
1841 //-----------------------------------------------------------------
1842 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1843 //-----------------------------------------------------------------
1849 void AliTPC::SetSide(Float_t side=0.)
1851 // choice of the TPC side
1856 //____________________________________________________________________________
1857 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
1858 Float_t p2,Float_t p3)
1861 // gax mixture definition
1875 //_____________________________________________________________________________
1877 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
1880 // electron transport taking into account:
1882 // 2.ExB at the wires
1883 // 3. nonisochronity
1885 // xyz and index must be already transformed to system 1
1888 fTPCParam->Transform1to2(xyz,index);
1891 Float_t driftl=xyz[2];
1892 if(driftl<0.01) driftl=0.01;
1893 driftl=TMath::Sqrt(driftl);
1894 Float_t sigT = driftl*(fTPCParam->GetDiffT());
1895 Float_t sigL = driftl*(fTPCParam->GetDiffL());
1896 xyz[0]=gRandom->Gaus(xyz[0],sigT);
1897 xyz[1]=gRandom->Gaus(xyz[1],sigT);
1898 xyz[2]=gRandom->Gaus(xyz[2],sigL);
1902 if (fTPCParam->GetMWPCReadout()==kTRUE){
1904 fTPCParam->Transform2to2NearestWire(xyz,index);
1905 Float_t dx=xyz[0]-x1;
1906 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
1908 //add nonisochronity (not implemented yet)
1912 ClassImp(AliTPCdigit)
1914 //_____________________________________________________________________________
1915 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
1919 // Creates a TPC digit object
1921 fSector = digits[0];
1922 fPadRow = digits[1];
1925 fSignal = digits[4];
1931 //_____________________________________________________________________________
1932 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
1936 // Creates a TPC hit object
1947 //________________________________________________________________________
1948 // Additional code because of the AliTPCTrackHits
1950 void AliTPC::MakeBranch2(Option_t *option,char *file)
1953 // Create a new branch in the current Root Tree
1954 // The branch of fHits is automatically split
1955 // MI change 14.09.2000
1956 if (fHitType&2==0) return;
1957 char branchname[10];
1958 sprintf(branchname,"%s2",GetName());
1960 // Get the pointer to the header
1961 char *cH = strstr(option,"H");
1963 if (fTrackHits && gAlice->TreeH() && cH) {
1964 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHits,
1965 gAlice->TreeH(),fBufferSize,1);
1966 gAlice->TreeH()->GetListOfBranches()->Add(branch);
1967 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
1969 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
1970 TDirectory *wd = gDirectory;
1972 TIter next( b->GetListOfBranches());
1973 while ((b=(TBranch*)next())) {
1977 cout << "Diverting branch " << branchname << " to file " << file << endl;
1982 void AliTPC::SetTreeAddress()
1984 if (fHitType&1) AliDetector::SetTreeAddress();
1985 if (fHitType&2) SetTreeAddress2();
1988 void AliTPC::SetTreeAddress2()
1991 // Set branch address for the TrackHits Tree
1994 char branchname[20];
1995 sprintf(branchname,"%s2",GetName());
1997 // Branch address for hit tree
1998 TTree *treeH = gAlice->TreeH();
2000 branch = treeH->GetBranch(branchname);
2001 if (branch) branch->SetAddress(&fTrackHits);
2005 void AliTPC::FinishPrimary()
2007 if (fTrackHits) fTrackHits->FlushHitStack();
2011 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2014 // add hit to the list
2017 int primary = gAlice->GetPrimary(track);
2018 gAlice->Particle(primary)->SetBit(kKeepBit);
2022 gAlice->FlagTrack(track);
2024 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2025 //if (hit->fTrack!=rtrack)
2026 // cout<<"bad track number\n";
2028 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2029 hits[1],hits[2],(Int_t)hits[3]);
2032 void AliTPC::ResetHits()
2034 if (fHitType&1) AliDetector::ResetHits();
2035 if (fHitType&2) ResetHits2();
2038 void AliTPC::ResetHits2()
2042 if (fTrackHits) fTrackHits->Clear();
2045 AliHit* AliTPC::FirstHit(Int_t track)
2047 if (fHitType&2) return FirstHit2(track);
2048 return AliDetector::FirstHit(track);
2050 AliHit* AliTPC::NextHit()
2052 if (fHitType&2) return NextHit2();
2053 return AliDetector::NextHit();
2056 AliHit* AliTPC::FirstHit2(Int_t track)
2059 // Initialise the hit iterator
2060 // Return the address of the first hit for track
2061 // If track>=0 the track is read from disk
2062 // while if track<0 the first hit of the current
2063 // track is returned
2066 gAlice->ResetHits();
2067 gAlice->TreeH()->GetEvent(track);
2071 fTrackHits->First();
2072 return fTrackHits->GetHit();
2077 AliHit* AliTPC::NextHit2()
2080 //Return the next hit for the current track
2084 return fTrackHits->GetHit();
2090 void AliTPC::LoadPoints(Int_t)
2094 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2097 if(fHitType==1) AliDetector::LoadPoints(a);
2098 else LoadPoints2(a);
2105 void AliTPC::RemapTrackHitIDs(Int_t *map)
2107 if (!fTrackHits) return;
2108 AliObjectArray * arr = fTrackHits->fTrackHitsInfo;
2109 for (UInt_t i=0;i<arr->GetSize();i++){
2110 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2111 info->fTrackID = map[info->fTrackID];
2117 //_____________________________________________________________________________
2118 void AliTPC::LoadPoints2(Int_t)
2121 // Store x, y, z of all hits in memory
2123 if (fTrackHits == 0) return;
2125 Int_t nhits = fTrackHits->GetEntriesFast();
2126 if (nhits == 0) return;
2127 Int_t tracks = gAlice->GetNtrack();
2128 if (fPoints == 0) fPoints = new TObjArray(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++) {
2150 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2151 trk=ahit->GetTrack();
2152 if(ntrk[trk]==limi[trk]) {
2154 // Initialise a new track
2155 fp=new Float_t[3*(limi[trk]+chunk)];
2157 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2158 delete [] coor[trk];
2165 fp[3*ntrk[trk] ] = ahit->X();
2166 fp[3*ntrk[trk]+1] = ahit->Y();
2167 fp[3*ntrk[trk]+2] = ahit->Z();
2172 for(trk=0; trk<tracks; ++trk) {
2174 points = new AliPoints();
2175 points->SetMarkerColor(GetMarkerColor());
2176 points->SetMarkerSize(GetMarkerSize());
2177 points->SetDetector(this);
2178 points->SetParticle(trk);
2179 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2180 fPoints->AddAt(points,trk);
2181 delete [] coor[trk];
2191 //_____________________________________________________________________________
2192 void AliTPC::LoadPoints3(Int_t)
2195 // Store x, y, z of all hits in memory
2196 // - only intersection point with pad row
2197 if (fTrackHits == 0) return;
2199 Int_t nhits = fTrackHits->GetEntriesFast();
2200 if (nhits == 0) return;
2201 Int_t tracks = gAlice->GetNtrack();
2202 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2203 fPoints->Expand(2*tracks);
2206 Int_t *ntrk=new Int_t[tracks];
2207 Int_t *limi=new Int_t[tracks];
2208 Float_t **coor=new Float_t*[tracks];
2209 for(Int_t i=0;i<tracks;i++) {
2215 AliPoints *points = 0;
2218 Int_t chunk=nhits/4+1;
2220 // Loop over all the hits and store their position
2222 ahit = FirstHit2(-1);
2223 //for (Int_t hit=0;hit<nhits;hit++) {
2227 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2228 trk=ahit->GetTrack();
2229 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2230 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2231 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2232 if (currentrow!=lastrow){
2233 lastrow = currentrow;
2234 //later calculate intersection point
2235 if(ntrk[trk]==limi[trk]) {
2237 // Initialise a new track
2238 fp=new Float_t[3*(limi[trk]+chunk)];
2240 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2241 delete [] coor[trk];
2248 fp[3*ntrk[trk] ] = ahit->X();
2249 fp[3*ntrk[trk]+1] = ahit->Y();
2250 fp[3*ntrk[trk]+2] = ahit->Z();
2257 for(trk=0; trk<tracks; ++trk) {
2259 points = new AliPoints();
2260 points->SetMarkerColor(GetMarkerColor()+1);
2261 points->SetMarkerStyle(5);
2262 points->SetMarkerSize(0.2);
2263 points->SetDetector(this);
2264 points->SetParticle(trk);
2265 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2266 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2267 fPoints->AddAt(points,tracks+trk);
2268 delete [] coor[trk];
2279 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2283 //fill clones array with intersection of current point with the
2288 const Int_t kcmaxhits=30000;
2289 TVector * xxxx = new TVector(kcmaxhits*4);
2290 TVector & xxx = *xxxx;
2291 Int_t maxhits = kcmaxhits;
2294 AliTPChit * tpcHit=0;
2295 tpcHit = (AliTPChit*)FirstHit2(-1);
2296 Int_t currentIndex=0;
2297 Int_t lastrow=-1; //last writen row
2300 if (tpcHit==0) continue;
2301 sector=tpcHit->fSector; // sector number
2302 ipart=tpcHit->Track();
2306 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2307 Int_t index[3]={1,sector,0};
2308 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2309 if (currentrow<0) continue;
2310 if (lastrow<0) lastrow=currentrow;
2311 if (currentrow==lastrow){
2312 if ( currentIndex>=maxhits){
2314 xxx.ResizeTo(4*maxhits);
2316 xxx(currentIndex*4)=x[0];
2317 xxx(currentIndex*4+1)=x[1];
2318 xxx(currentIndex*4+2)=x[2];
2319 xxx(currentIndex*4+3)=tpcHit->fQ;
2323 if (currentIndex>2){
2335 for (Int_t index=0;index<currentIndex;index++){
2337 x=x2=x3=x4=xxx(index*4);
2345 sumy+=xxx(index*4+1);
2346 sumxy+=xxx(index*4+1)*x;
2347 sumx2y+=xxx(index*4+1)*x2;
2348 sumz+=xxx(index*4+2);
2349 sumxz+=xxx(index*4+2)*x;
2350 sumx2z+=xxx(index*4+2)*x2;
2351 sumq+=xxx(index*4+3);
2353 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2354 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2355 sumx2*(sumx*sumx3-sumx2*sumx2);
2357 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2358 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2359 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2360 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2362 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2363 sumx2*(sumx*sumx2y-sumx2*sumxy);
2364 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2365 sumx2*(sumx*sumx2z-sumx2*sumxz);
2367 Float_t y=detay/det+centralPad;
2368 Float_t z=detaz/det;
2369 Float_t by=detby/det; //y angle
2370 Float_t bz=detbz/det; //z angle
2371 sumy/=Float_t(currentIndex);
2372 sumz/=Float_t(currentIndex);
2374 AliComplexCluster cl;
2380 cl.fTracks[0]=ipart;
2382 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2383 if (row!=0) row->InsertCluster(&cl);
2386 } //end of calculating cluster for given row
2388 } // end of loop over hits