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.31 2001/03/12 08:21:50 kowal2
19 Corrected C++ bug in the material definitions
21 Revision 1.30 2001/03/01 17:34:47 kowal2
22 Correction due to the accuracy problem
24 Revision 1.29 2001/02/28 16:34:40 kowal2
25 Protection against nonphysical values of the avalanche size,
28 Revision 1.28 2001/01/26 19:57:19 hristov
29 Major upgrade of AliRoot code
31 Revision 1.27 2001/01/13 17:29:33 kowal2
32 Sun compiler correction
34 Revision 1.26 2001/01/10 07:59:43 kowal2
35 Corrections to load points from the noncompressed hits.
37 Revision 1.25 2000/11/02 07:25:31 kowal2
38 Changes due to the new hit structure.
41 Revision 1.24 2000/10/05 16:06:09 kowal2
42 Forward declarations. Changes due to a new class AliComplexCluster.
44 Revision 1.23 2000/10/02 21:28:18 fca
45 Removal of useless dependecies via forward declarations
47 Revision 1.22 2000/07/10 20:57:39 hristov
48 Update of TPC code and macros by M.Kowalski
50 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
51 Changes to obey the coding rules
53 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
54 Splitted from AliTPCtracking
56 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
57 Changed parameter settings
59 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
61 Defaults loaded automatically (hard-wired)
62 Optional parameters can be set via macro called in the constructor
64 Revision 1.19 2000/04/18 19:00:59 fca
65 Small bug fixes to TPC files
67 Revision 1.18 2000/04/17 09:37:33 kowal2
68 removed obsolete AliTPCDigitsDisplay.C
70 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
72 New, experimental data structure from M. Ivanov
73 New tracking algorithm
74 Different pad geometry for different sectors
75 Digitization rewritten
77 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
78 Not used anymore - removed
80 Revision 1.17 2000/01/19 17:17:30 fca
81 Introducing a list of lists of hits -- more hits allowed for detector now
83 Revision 1.16 1999/11/05 09:29:23 fca
84 Accept only signals > 0
86 Revision 1.15 1999/10/08 06:26:53 fca
87 Removed ClustersIndex - not used anymore
89 Revision 1.14 1999/09/29 09:24:33 fca
90 Introduction of the Copyright and cvs Log
94 ///////////////////////////////////////////////////////////////////////////////
96 // Time Projection Chamber //
97 // This class contains the basic functions for the Time Projection Chamber //
98 // detector. Functions specific to one particular geometry are //
99 // contained in the derived classes //
103 <img src="picts/AliTPCClass.gif">
108 ///////////////////////////////////////////////////////////////////////////////
116 #include <TGeometry.h>
119 #include <TObjectTable.h>
120 #include "TParticle.h"
124 #include <iostream.h>
131 #include "AliTPCParamSR.h"
132 #include "AliTPCPRF2D.h"
133 #include "AliTPCRF1D.h"
134 #include "AliDigits.h"
135 #include "AliSimDigits.h"
136 #include "AliTPCTrackHits.h"
137 #include "AliPoints.h"
138 #include "AliArrayBranch.h"
141 #include "AliTPCDigitsArray.h"
142 #include "AliComplexCluster.h"
143 #include "AliClusters.h"
144 #include "AliTPCClustersRow.h"
145 #include "AliTPCClustersArray.h"
147 #include "AliTPCcluster.h"
148 #include "AliTPCclusterer.h"
149 #include "AliTPCtracker.h"
151 #include <TInterpreter.h>
158 //_____________________________________________________________________________
162 // Default constructor
177 //_____________________________________________________________________________
178 AliTPC::AliTPC(const char *name, const char *title)
179 : AliDetector(name,title)
182 // Standard constructor
186 // Initialise arrays of hits and digits
187 fHits = new TClonesArray("AliTPChit", 176);
188 gAlice->AddHitList(fHits);
193 fTrackHits = new AliTPCTrackHits; //MI - 13.09.2000
194 fTrackHits->SetHitPrecision(0.002);
195 fTrackHits->SetStepPrecision(0.003);
196 fTrackHits->SetMaxDistance(100);
200 // Initialise counters
206 // Initialise color attributes
207 SetMarkerColor(kYellow);
210 // Set TPC parameters
213 if (!strcmp(title,"Default")) {
214 fTPCParam = new AliTPCParamSR;
216 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
222 //_____________________________________________________________________________
233 delete fTrackHits; //MI 15.09.2000
236 //_____________________________________________________________________________
237 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
240 // Add a hit to the list
242 // TClonesArray &lhits = *fHits;
243 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
245 TClonesArray &lhits = *fHits;
246 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
249 AddHit2(track,vol,hits);
252 //_____________________________________________________________________________
253 void AliTPC::BuildGeometry()
257 // Build TPC ROOT TNode geometry for the event display
262 const int kColorTPC=19;
263 char name[5], title[25];
264 const Double_t kDegrad=TMath::Pi()/180;
265 const Double_t kRaddeg=180./TMath::Pi();
268 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
269 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
271 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
272 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
274 Int_t nLo = fTPCParam->GetNInnerSector()/2;
275 Int_t nHi = fTPCParam->GetNOuterSector()/2;
277 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
278 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
279 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
280 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
283 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
284 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
290 // Get ALICE top node
293 nTop=gAlice->GetGeometry()->GetNode("alice");
297 rl = fTPCParam->GetInnerRadiusLow();
298 ru = fTPCParam->GetInnerRadiusUp();
302 sprintf(name,"LS%2.2d",i);
304 sprintf(title,"TPC low sector %3d",i);
307 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
308 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
309 tubs->SetNumberOfDivisions(1);
311 nNode = new TNode(name,title,name,0,0,0,"");
312 nNode->SetLineColor(kColorTPC);
318 rl = fTPCParam->GetOuterRadiusLow();
319 ru = fTPCParam->GetOuterRadiusUp();
322 sprintf(name,"US%2.2d",i);
324 sprintf(title,"TPC upper sector %d",i);
326 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
327 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
328 tubs->SetNumberOfDivisions(1);
330 nNode = new TNode(name,title,name,0,0,0,"");
331 nNode->SetLineColor(kColorTPC);
337 //_____________________________________________________________________________
338 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
341 // Calculate distance from TPC to mouse on the display
347 void AliTPC::Clusters2Tracks(TFile *of) {
348 //-----------------------------------------------------------------
349 // This is a track finder.
350 //-----------------------------------------------------------------
351 AliTPCtracker::Clusters2Tracks(fTPCParam,of);
354 //_____________________________________________________________________________
355 void AliTPC::CreateMaterials()
357 //-----------------------------------------------
358 // Create Materials for for TPC simulations
359 //-----------------------------------------------
361 //-----------------------------------------------------------------
362 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
363 //-----------------------------------------------------------------
365 Int_t iSXFLD=gAlice->Field()->Integ();
366 Float_t sXMGMX=gAlice->Field()->Max();
368 Float_t amat[5]; // atomic numbers
369 Float_t zmat[5]; // z
370 Float_t wmat[5]; // proportions
376 //***************** Gases *************************
378 //-------------------------------------------------
380 //-------------------------------------------------
391 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
401 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
404 //--------------------------------------------------------------
406 //--------------------------------------------------------------
423 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
425 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
440 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
442 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
458 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
460 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
462 //----------------------------------------------------------------
463 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
464 //----------------------------------------------------------------
470 Float_t rho,absl,X0,buf[1];
474 for(nc = 0;nc<fNoComp;nc++)
477 // retrive material constants
479 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
484 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
486 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
487 density += fMixtProp[nc]*rho; // density of the mixture
491 // mixture proportions by weight!
493 for(nc = 0;nc<fNoComp;nc++)
496 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
498 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
499 apure[nnc] : amol[nnc])/am;
503 // Drift gases 1 - nonsensitive, 2 - sensitive
505 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
506 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
515 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
518 //----------------------------------------------------------------------
520 //----------------------------------------------------------------------
542 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
564 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
582 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
600 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
602 // G10 60% SiO2 + 40% epoxy, I use A and Z for SiO2
615 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz
617 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
619 AliMaterial(39,"G10",amat[0],zmat[0],density,999.,999.);
628 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
637 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
646 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
664 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
683 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
687 //----------------------------------------------------------
688 // tracking media for gases
689 //----------------------------------------------------------
691 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
692 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
693 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
694 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
696 //-----------------------------------------------------------
697 // tracking media for solids
698 //-----------------------------------------------------------
700 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
701 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
702 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
703 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
704 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
705 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
706 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
707 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
708 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
709 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
714 void AliTPC::Digits2Clusters(TFile *of)
716 //-----------------------------------------------------------------
717 // This is a simple cluster finder.
718 //-----------------------------------------------------------------
719 AliTPCclusterer::Digits2Clusters(fTPCParam,of);
722 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
723 extern Double_t SigmaZ2(Double_t, Double_t);
724 //_____________________________________________________________________________
725 void AliTPC::Hits2Clusters(TFile *of)
727 //--------------------------------------------------------
728 // TPC simple cluster generator from hits
729 // obtained from the TPC Fast Simulator
730 // The point errors are taken from the parametrization
731 //--------------------------------------------------------
733 //-----------------------------------------------------------------
734 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
735 //-----------------------------------------------------------------
736 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
737 // Jouri.Belikov@cern.ch
738 //----------------------------------------------------------------
740 /////////////////////////////////////////////////////////////////////////////
742 //---------------------------------------------------------------------
743 // ALICE TPC Cluster Parameters
744 //--------------------------------------------------------------------
748 // Cluster width in rphi
749 const Float_t kACrphi=0.18322;
750 const Float_t kBCrphi=0.59551e-3;
751 const Float_t kCCrphi=0.60952e-1;
752 // Cluster width in z
753 const Float_t kACz=0.19081;
754 const Float_t kBCz=0.55938e-3;
755 const Float_t kCCz=0.30428;
757 TDirectory *savedir=gDirectory;
760 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
765 printf("AliTPCParam MUST be created firstly\n");
769 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
771 TParticle *particle; // pointer to a given particle
772 AliTPChit *tpcHit; // pointer to a sigle TPC hit
776 Float_t pl,pt,tanth,rpad,ratio;
779 //---------------------------------------------------------------
780 // Get the access to the tracks
781 //---------------------------------------------------------------
783 TTree *tH = gAlice->TreeH();
784 Stat_t ntracks = tH->GetEntries();
786 //Switch to the output file
789 fTPCParam->Write(fTPCParam->GetTitle());
790 AliTPCClustersArray carray;
791 carray.Setup(fTPCParam);
792 carray.SetClusterType("AliTPCcluster");
795 Int_t nclusters=0; //cluster counter
797 //------------------------------------------------------------
798 // Loop over all sectors (72 sectors for 20 deg
799 // segmentation for both lower and upper sectors)
800 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
801 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
803 // First cluster for sector 0 starts at "0"
804 //------------------------------------------------------------
806 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
808 fTPCParam->AdjustCosSin(isec,cph,sph);
810 //------------------------------------------------------------
812 //------------------------------------------------------------
814 for(Int_t track=0;track<ntracks;track++){
818 // Get number of the TPC hits
820 // nhits=fHits->GetEntriesFast();
823 tpcHit = (AliTPChit*)FirstHit(-1);
827 // for(Int_t hit=0;hit<nhits;hit++){
828 //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
832 if (tpcHit->fQ == 0.) {
833 tpcHit = (AliTPChit*) NextHit();
834 continue; //information about track (I.Belikov)
836 sector=tpcHit->fSector; // sector number
839 // if(sector != isec) continue; //terminate iteration
842 tpcHit = (AliTPChit*) NextHit();
845 ipart=tpcHit->Track();
846 particle=gAlice->Particle(ipart);
849 if(pt < 1.e-9) pt=1.e-9;
851 tanth = TMath::Abs(tanth);
852 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
853 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
855 // space-point resolutions
857 sigmaRphi=SigmaY2(rpad,tanth,pt);
858 sigmaZ =SigmaZ2(rpad,tanth );
862 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
863 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
865 // temporary protection
867 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
868 if(sigmaZ < 0.) sigmaZ=0.4e-3;
869 if(clRphi < 0.) clRphi=2.5e-3;
870 if(clZ < 0.) clZ=2.5e-5;
875 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
876 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
878 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
879 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
880 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
881 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
882 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
883 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
884 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
885 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
886 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
887 xyz[2]=sigmaRphi; // fSigmaY2
888 xyz[3]=sigmaZ; // fSigmaZ2
889 xyz[4]=tpcHit->fQ; // q
891 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
892 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
894 Int_t tracks[3]={tpcHit->Track(), -1, -1};
895 AliTPCcluster cluster(tracks,xyz);
897 clrow->InsertCluster(&cluster); nclusters++;
899 tpcHit = (AliTPChit*)NextHit();
902 } // end of loop over hits
904 } // end of loop over tracks
906 Int_t nrows=fTPCParam->GetNRow(isec);
907 for (Int_t irow=0; irow<nrows; irow++) {
908 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
909 if (!clrow) continue;
910 carray.StoreRow(isec,irow);
911 carray.ClearRow(isec,irow);
914 } // end of loop over sectors
916 cerr<<"Number of made clusters : "<<nclusters<<" \n";
918 carray.GetTree()->Write();
920 savedir->cd(); //switch back to the input file
924 //_________________________________________________________________
925 void AliTPC::Hits2ExactClustersSector(Int_t isec)
927 //--------------------------------------------------------
928 //calculate exact cross point of track and given pad row
929 //resulting values are expressed in "digit" coordinata
930 //--------------------------------------------------------
932 //-----------------------------------------------------------------
933 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
934 //-----------------------------------------------------------------
936 if (fClustersArray==0){
940 TParticle *particle; // pointer to a given particle
941 AliTPChit *tpcHit; // pointer to a sigle TPC hit
944 const Int_t kcmaxhits=30000;
945 TVector * xxxx = new TVector(kcmaxhits*4);
946 TVector & xxx = *xxxx;
947 Int_t maxhits = kcmaxhits;
948 //construct array for each padrow
949 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
950 fClustersArray->CreateRow(isec,i);
952 //---------------------------------------------------------------
953 // Get the access to the tracks
954 //---------------------------------------------------------------
956 TTree *tH = gAlice->TreeH();
957 Stat_t ntracks = tH->GetEntries();
958 Int_t npart = gAlice->GetNtrack();
960 //------------------------------------------------------------
962 //------------------------------------------------------------
964 for(Int_t track=0;track<ntracks;track++){
968 // Get number of the TPC hits and a pointer
971 nhits=fHits->GetEntriesFast();
975 Int_t currentIndex=0;
976 Int_t lastrow=-1; //last writen row
977 for(Int_t hit=0;hit<nhits;hit++){
978 tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
979 if (tpcHit==0) continue;
980 sector=tpcHit->fSector; // sector number
981 if(sector != isec) continue;
982 ipart=tpcHit->Track();
983 if (ipart<npart) particle=gAlice->Particle(ipart);
987 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
988 Int_t index[3]={1,isec,0};
989 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
990 if (currentrow<0) continue;
991 if (lastrow<0) lastrow=currentrow;
992 if (currentrow==lastrow){
993 if ( currentIndex>=maxhits){
995 xxx.ResizeTo(4*maxhits);
997 xxx(currentIndex*4)=x[0];
998 xxx(currentIndex*4+1)=x[1];
999 xxx(currentIndex*4+2)=x[2];
1000 xxx(currentIndex*4+3)=tpcHit->fQ;
1004 if (currentIndex>2){
1016 for (Int_t index=0;index<currentIndex;index++){
1018 x=x2=x3=x4=xxx(index*4);
1026 sumy+=xxx(index*4+1);
1027 sumxy+=xxx(index*4+1)*x;
1028 sumx2y+=xxx(index*4+1)*x2;
1029 sumz+=xxx(index*4+2);
1030 sumxz+=xxx(index*4+2)*x;
1031 sumx2z+=xxx(index*4+2)*x2;
1032 sumq+=xxx(index*4+3);
1034 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1035 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1036 sumx2*(sumx*sumx3-sumx2*sumx2);
1038 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1039 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1040 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1041 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1043 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1044 sumx2*(sumx*sumx2y-sumx2*sumxy);
1045 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1046 sumx2*(sumx*sumx2z-sumx2*sumxz);
1048 Float_t y=detay/det+centralPad;
1049 Float_t z=detaz/det;
1050 Float_t by=detby/det; //y angle
1051 Float_t bz=detbz/det; //z angle
1052 sumy/=Float_t(currentIndex);
1053 sumz/=Float_t(currentIndex);
1054 AliComplexCluster cl;
1060 cl.fTracks[0]=ipart;
1062 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1063 if (row!=0) row->InsertCluster(&cl);
1066 } //end of calculating cluster for given row
1070 } // end of loop over hits
1071 } // end of loop over tracks
1072 //write padrows to tree
1073 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1074 fClustersArray->StoreRow(isec,ii);
1075 fClustersArray->ClearRow(isec,ii);
1080 //___________________________________________
1081 void AliTPC::SDigits2Digits()
1083 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1084 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1085 AliTPCPRF2D * prfouter = new AliTPCPRF2D;
1086 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1088 TDirectory *cwd = gDirectory;
1089 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1090 rf->SetOffset(3*param->GetZSigma());
1092 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1094 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n";
1097 prfinner->Read("prf_07504_Gati_056068_d02");
1098 prfouter->Read("prf_10006_Gati_047051_d03");
1102 param->SetInnerPRF(prfinner);
1103 param->SetOuterPRF(prfouter);
1104 param->SetTimeRF(rf);
1108 cerr<<"Digitizing TPC...\n";
1110 //setup TPCDigitsArray
1111 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1112 arr->SetClass("AliSimDigits");
1116 arr->MakeTree(fDigitsFile);
1118 SetDigitsArray(arr);
1122 // Hits2DigitsSector(1);
1123 // Hits2DigitsSector(2);
1124 // Hits2DigitsSector(3);
1125 // Hits2DigitsSector(1+18);
1126 // Hits2DigitsSector(2+18);
1127 // Hits2DigitsSector(3+18);
1129 // Hits2DigitsSector(36+1);
1130 // Hits2DigitsSector(36+2);
1131 // Hits2DigitsSector(36+3);
1132 // Hits2DigitsSector(36+1+18);
1133 // Hits2DigitsSector(36+2+18);
1134 // Hits2DigitsSector(36+3+18);
1139 sprintf(treeName,"TreeD_%s",param->GetTitle());
1140 GetDigitsArray()->GetTree()->Write(treeName,TObject::kOverwrite);
1143 //__________________________________________________________________
1144 void AliTPC::Hits2Digits()
1146 //----------------------------------------------------
1147 // Loop over all sectors
1148 //----------------------------------------------------
1151 printf("AliTPCParam MUST be created firstly\n");
1155 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1160 //_____________________________________________________________________________
1161 void AliTPC::Hits2DigitsSector(Int_t isec)
1163 //-------------------------------------------------------------------
1164 // TPC conversion from hits to digits.
1165 //-------------------------------------------------------------------
1167 //-----------------------------------------------------------------
1168 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1169 //-----------------------------------------------------------------
1171 //-------------------------------------------------------
1172 // Get the access to the track hits
1173 //-------------------------------------------------------
1176 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1177 Stat_t ntracks = tH->GetEntries();
1181 //-------------------------------------------
1182 // Only if there are any tracks...
1183 //-------------------------------------------
1187 printf("*** Processing sector number %d ***\n",isec);
1189 Int_t nrows =fTPCParam->GetNRow(isec);
1191 row= new TObjArray* [nrows];
1193 MakeSector(isec,nrows,tH,ntracks,row);
1195 //--------------------------------------------------------
1196 // Digitize this sector, row by row
1197 // row[i] is the pointer to the TObjArray of TVectors,
1198 // each one containing electrons accepted on this
1199 // row, assigned into tracks
1200 //--------------------------------------------------------
1204 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1206 for (i=0;i<nrows;i++){
1208 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1210 DigitizeRow(i,isec,row);
1212 fDigitsArray->StoreRow(isec,i);
1214 Int_t ndig = dig->GetDigitSize();
1216 printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1218 fDigitsArray->ClearRow(isec,i);
1221 } // end of the sector digitization
1223 for(i=0;i<nrows;i++){
1228 delete [] row; // delete the array of pointers to TObjArray-s
1232 } // end of Hits2DigitsSector
1235 //_____________________________________________________________________________
1236 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1238 //-----------------------------------------------------------
1239 // Single row digitization, coupling from the neighbouring
1240 // rows taken into account
1241 //-----------------------------------------------------------
1243 //-----------------------------------------------------------------
1244 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1245 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1246 //-----------------------------------------------------------------
1249 Float_t zerosup = fTPCParam->GetZeroSup();
1250 Int_t nrows =fTPCParam->GetNRow(isec);
1251 fCurrentIndex[1]= isec;
1254 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1255 Int_t nofTbins = fTPCParam->GetMaxTBin();
1256 Int_t indexRange[4];
1258 // Integrated signal for this row
1259 // and a single track signal
1261 TMatrix *m1 = new TMatrix(0,nofPads,0,nofTbins); // integrated
1262 TMatrix *m2 = new TMatrix(0,nofPads,0,nofTbins); // single
1264 TMatrix &total = *m1;
1266 // Array of pointers to the label-signal list
1268 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1269 Float_t **pList = new Float_t* [nofDigits];
1273 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1277 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1278 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1279 for (Int_t row= row1;row<=row2;row++){
1280 Int_t nTracks= rows[row]->GetEntries();
1281 for (i1=0;i1<nTracks;i1++){
1282 fCurrentIndex[2]= row;
1283 fCurrentIndex[3]=irow;
1285 m2->Zero(); // clear single track signal matrix
1286 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1287 GetList(trackLabel,nofPads,m2,indexRange,pList);
1289 else GetSignal(rows[row],i1,0,m1,indexRange);
1295 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1296 for(Int_t ip=0;ip<nofPads;ip++){
1297 for(Int_t it=0;it<nofTbins;it++){
1299 Float_t q = total(ip,it);
1301 Int_t gi =it*nofPads+ip; // global index
1303 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
1307 if(q <=zerosup) continue; // do not fill zeros
1308 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1311 // "real" signal or electronic noise (list = -1)?
1314 for(Int_t j1=0;j1<3;j1++){
1315 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1320 <A NAME="AliDigits"></A>
1321 using of AliDigits object
1324 dig->SetDigitFast((Short_t)q,it,ip);
1325 if (fDigitsArray->IsSimulated())
1327 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1328 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1329 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1333 } // end of loop over time buckets
1334 } // end of lop over pads
1337 // This row has been digitized, delete nonused stuff
1340 for(lp=0;lp<nofDigits;lp++){
1341 if(pList[lp]) delete [] pList[lp];
1350 } // end of DigitizeRow
1352 //_____________________________________________________________________________
1354 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, TMatrix *m1, TMatrix *m2,
1358 //---------------------------------------------------------------
1359 // Calculates 2-D signal (pad,time) for a single track,
1360 // returns a pointer to the signal matrix and the track label
1361 // No digitization is performed at this level!!!
1362 //---------------------------------------------------------------
1364 //-----------------------------------------------------------------
1365 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1366 // Modified: Marian Ivanov
1367 //-----------------------------------------------------------------
1371 tv = (TVector*)p1->At(ntr); // pointer to a track
1374 Float_t label = v(0);
1375 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1377 Int_t nElectrons = (tv->GetNrows()-1)/4;
1378 indexRange[0]=9999; // min pad
1379 indexRange[1]=-1; // max pad
1380 indexRange[2]=9999; //min time
1381 indexRange[3]=-1; // max time
1383 // Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1385 TMatrix &signal = *m1;
1386 TMatrix &total = *m2;
1388 // Loop over all electrons
1390 for(Int_t nel=0; nel<nElectrons; nel++){
1392 Float_t aval = v(idx+4);
1393 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1394 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1395 Int_t n = fTPCParam->CalcResponse(xyz,fCurrentIndex,fCurrentIndex[3]);
1397 if (n>0) for (Int_t i =0; i<n; i++){
1398 Int_t *index = fTPCParam->GetResBin(i);
1399 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1400 if ( ( pad<(fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]))) && (pad>0)) {
1401 Int_t time=index[2];
1402 Float_t weight = fTPCParam->GetResWeight(i); //we normalise response to ADC channel
1403 weight *= eltoadcfac;
1405 if (m1!=0) signal(pad,time)+=weight;
1406 total(pad,time)+=weight;
1407 indexRange[0]=TMath::Min(indexRange[0],pad);
1408 indexRange[1]=TMath::Max(indexRange[1],pad);
1409 indexRange[2]=TMath::Min(indexRange[2],time);
1410 indexRange[3]=TMath::Max(indexRange[3],time);
1413 } // end of loop over electrons
1415 return label; // returns track label when finished
1418 //_____________________________________________________________________________
1419 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *indexRange,
1422 //----------------------------------------------------------------------
1423 // Updates the list of tracks contributing to digits for a given row
1424 //----------------------------------------------------------------------
1426 //-----------------------------------------------------------------
1427 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1428 //-----------------------------------------------------------------
1430 TMatrix &signal = *m;
1432 // lop over nonzero digits
1434 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1435 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1438 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1440 if(signal(ip,it)<0.5) continue;
1443 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1445 if(!pList[globalIndex]){
1448 // Create new list (6 elements - 3 signals and 3 labels),
1451 pList[globalIndex] = new Float_t [6];
1455 *pList[globalIndex] = -1.;
1456 *(pList[globalIndex]+1) = -1.;
1457 *(pList[globalIndex]+2) = -1.;
1458 *(pList[globalIndex]+3) = -1.;
1459 *(pList[globalIndex]+4) = -1.;
1460 *(pList[globalIndex]+5) = -1.;
1463 *pList[globalIndex] = label;
1464 *(pList[globalIndex]+3) = signal(ip,it);
1468 // check the signal magnitude
1470 Float_t highest = *(pList[globalIndex]+3);
1471 Float_t middle = *(pList[globalIndex]+4);
1472 Float_t lowest = *(pList[globalIndex]+5);
1475 // compare the new signal with already existing list
1478 if(signal(ip,it)<lowest) continue; // neglect this track
1482 if (signal(ip,it)>highest){
1483 *(pList[globalIndex]+5) = middle;
1484 *(pList[globalIndex]+4) = highest;
1485 *(pList[globalIndex]+3) = signal(ip,it);
1487 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1488 *(pList[globalIndex]+1) = *pList[globalIndex];
1489 *pList[globalIndex] = label;
1491 else if (signal(ip,it)>middle){
1492 *(pList[globalIndex]+5) = middle;
1493 *(pList[globalIndex]+4) = signal(ip,it);
1495 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1496 *(pList[globalIndex]+1) = label;
1499 *(pList[globalIndex]+5) = signal(ip,it);
1500 *(pList[globalIndex]+2) = label;
1504 } // end of loop over pads
1505 } // end of loop over time bins
1510 //___________________________________________________________________
1511 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1512 Stat_t ntracks,TObjArray **row)
1515 //-----------------------------------------------------------------
1516 // Prepares the sector digitization, creates the vectors of
1517 // tracks for each row of this sector. The track vector
1518 // contains the track label and the position of electrons.
1519 //-----------------------------------------------------------------
1521 //-----------------------------------------------------------------
1522 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1523 //-----------------------------------------------------------------
1525 Float_t gasgain = fTPCParam->GetGasGain();
1529 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1532 if (fHitType&2) branch = TH->GetBranch("TPC2");
1533 else branch = TH->GetBranch("TPC");
1536 //----------------------------------------------
1537 // Create TObjArray-s, one for each row,
1538 // each TObjArray will store the TVectors
1539 // of electrons, one TVector per each track.
1540 //----------------------------------------------
1542 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
1543 TVector **tracks = new TVector* [nrows]; //pointers to the track vectors
1544 for(i=0; i<nrows; i++){
1545 row[i] = new TObjArray;
1552 //--------------------------------------------------------------------
1553 // Loop over tracks, the "track" contains the full history
1554 //--------------------------------------------------------------------
1556 Int_t previousTrack,currentTrack;
1557 previousTrack = -1; // nothing to store so far!
1559 for(Int_t track=0;track<ntracks;track++){
1560 Bool_t isInSector=kTRUE;
1565 TBranch * br = TH->GetBranch("fTrackHitsInfo");
1566 br->GetEvent(track);
1567 AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
1568 for (UInt_t j=0;j<ar->GetSize();j++){
1569 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
1572 if (!isInSector) continue;
1574 branch->GetEntry(track); // get next track
1578 tpcHit = (AliTPChit*)FirstHit(-1);
1580 //--------------------------------------------------------------
1582 //--------------------------------------------------------------
1587 Int_t sector=tpcHit->fSector; // sector number
1588 // if(sector != isec) continue;
1590 tpcHit = (AliTPChit*) NextHit();
1594 currentTrack = tpcHit->Track(); // track number
1597 if(currentTrack != previousTrack){
1599 // store already filled fTrack
1601 for(i=0;i<nrows;i++){
1602 if(previousTrack != -1){
1603 if(nofElectrons[i]>0){
1604 TVector &v = *tracks[i];
1605 v(0) = previousTrack;
1606 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1607 row[i]->Add(tracks[i]);
1610 delete tracks[i]; // delete empty TVector
1616 tracks[i] = new TVector(481); // TVectors for the next fTrack
1618 } // end of loop over rows
1620 previousTrack=currentTrack; // update track label
1623 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1625 //---------------------------------------------------
1626 // Calculate the electron attachment probability
1627 //---------------------------------------------------
1630 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
1631 /fTPCParam->GetDriftV();
1633 Float_t attProb = fTPCParam->GetAttCoef()*
1634 fTPCParam->GetOxyCont()*time; // fraction!
1636 //-----------------------------------------------
1637 // Loop over electrons
1638 //-----------------------------------------------
1641 for(Int_t nel=0;nel<qI;nel++){
1642 // skip if electron lost due to the attachment
1643 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1648 // protection for the nonphysical avalanche size (10**6 maximum)
1650 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
1651 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
1654 TransportElectron(xyz,index); //MI change -august
1656 fTPCParam->GetPadRow(xyz,index); //MI change august
1657 rowNumber = index[2];
1658 //transform position to local digit coordinates
1659 //relative to nearest pad row
1660 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
1661 nofElectrons[rowNumber]++;
1662 //----------------------------------
1663 // Expand vector if necessary
1664 //----------------------------------
1665 if(nofElectrons[rowNumber]>120){
1666 Int_t range = tracks[rowNumber]->GetNrows();
1667 if((nofElectrons[rowNumber])>(range-1)/4){
1669 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
1673 TVector &v = *tracks[rowNumber];
1674 Int_t idx = 4*nofElectrons[rowNumber]-3;
1676 v(idx)= xyz[0]; // X - pad row coordinate
1677 v(idx+1)=xyz[1]; // Y - pad coordinate (along the pad-row)
1678 v(idx+2)=xyz[2]; // Z - time bin coordinate
1679 v(idx+3)=xyz[3]; // avalanche size
1680 } // end of loop over electrons
1682 tpcHit = (AliTPChit*)NextHit();
1684 } // end of loop over hits
1685 } // end of loop over tracks
1688 // store remaining track (the last one) if not empty
1691 for(i=0;i<nrows;i++){
1692 if(nofElectrons[i]>0){
1693 TVector &v = *tracks[i];
1694 v(0) = previousTrack;
1695 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1696 row[i]->Add(tracks[i]);
1705 delete [] nofElectrons;
1708 } // end of MakeSector
1711 //_____________________________________________________________________________
1715 // Initialise TPC detector after definition of geometry
1720 for(i=0;i<35;i++) printf("*");
1721 printf(" TPC_INIT ");
1722 for(i=0;i<35;i++) printf("*");
1725 for(i=0;i<80;i++) printf("*");
1729 //_____________________________________________________________________________
1730 void AliTPC::MakeBranch(Option_t* option, char *file)
1733 // Create Tree branches for the TPC.
1735 Int_t buffersize = 4000;
1736 char branchname[10];
1737 sprintf(branchname,"%s",GetName());
1739 AliDetector::MakeBranch(option,file);
1741 const char *d = strstr(option,"D");
1743 if (fDigits && gAlice->TreeD() && d) {
1744 gAlice->MakeBranchInTree(gAlice->TreeD(),
1745 branchname, &fDigits, buffersize, file) ;
1748 if (fHitType&2) MakeBranch2(option,file); // MI change 14.09.2000
1751 //_____________________________________________________________________________
1752 void AliTPC::ResetDigits()
1755 // Reset number of digits and the digits array for this detector
1758 if (fDigits) fDigits->Clear();
1761 //_____________________________________________________________________________
1762 void AliTPC::SetSecAL(Int_t sec)
1764 //---------------------------------------------------
1765 // Activate/deactivate selection for lower sectors
1766 //---------------------------------------------------
1768 //-----------------------------------------------------------------
1769 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1770 //-----------------------------------------------------------------
1775 //_____________________________________________________________________________
1776 void AliTPC::SetSecAU(Int_t sec)
1778 //----------------------------------------------------
1779 // Activate/deactivate selection for upper sectors
1780 //---------------------------------------------------
1782 //-----------------------------------------------------------------
1783 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1784 //-----------------------------------------------------------------
1789 //_____________________________________________________________________________
1790 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
1792 //----------------------------------------
1793 // Select active lower sectors
1794 //----------------------------------------
1796 //-----------------------------------------------------------------
1797 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1798 //-----------------------------------------------------------------
1808 //_____________________________________________________________________________
1809 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
1810 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
1811 Int_t s11 , Int_t s12)
1813 //--------------------------------
1814 // Select active upper sectors
1815 //--------------------------------
1817 //-----------------------------------------------------------------
1818 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1819 //-----------------------------------------------------------------
1835 //_____________________________________________________________________________
1836 void AliTPC::SetSens(Int_t sens)
1839 //-------------------------------------------------------------
1840 // Activates/deactivates the sensitive strips at the center of
1841 // the pad row -- this is for the space-point resolution calculations
1842 //-------------------------------------------------------------
1844 //-----------------------------------------------------------------
1845 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1846 //-----------------------------------------------------------------
1852 void AliTPC::SetSide(Float_t side=0.)
1854 // choice of the TPC side
1859 //____________________________________________________________________________
1860 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
1861 Float_t p2,Float_t p3)
1864 // gax mixture definition
1878 //_____________________________________________________________________________
1880 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
1883 // electron transport taking into account:
1885 // 2.ExB at the wires
1886 // 3. nonisochronity
1888 // xyz and index must be already transformed to system 1
1891 fTPCParam->Transform1to2(xyz,index);
1894 Float_t driftl=xyz[2];
1895 if(driftl<0.01) driftl=0.01;
1896 driftl=TMath::Sqrt(driftl);
1897 Float_t sigT = driftl*(fTPCParam->GetDiffT());
1898 Float_t sigL = driftl*(fTPCParam->GetDiffL());
1899 xyz[0]=gRandom->Gaus(xyz[0],sigT);
1900 xyz[1]=gRandom->Gaus(xyz[1],sigT);
1901 xyz[2]=gRandom->Gaus(xyz[2],sigL);
1905 if (fTPCParam->GetMWPCReadout()==kTRUE){
1907 fTPCParam->Transform2to2NearestWire(xyz,index);
1908 Float_t dx=xyz[0]-x1;
1909 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
1911 //add nonisochronity (not implemented yet)
1915 ClassImp(AliTPCdigit)
1917 //_____________________________________________________________________________
1918 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
1922 // Creates a TPC digit object
1924 fSector = digits[0];
1925 fPadRow = digits[1];
1928 fSignal = digits[4];
1934 //_____________________________________________________________________________
1935 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
1939 // Creates a TPC hit object
1950 //________________________________________________________________________
1951 // Additional code because of the AliTPCTrackHits
1953 void AliTPC::MakeBranch2(Option_t *option,char *file)
1956 // Create a new branch in the current Root Tree
1957 // The branch of fHits is automatically split
1958 // MI change 14.09.2000
1959 if (fHitType&2==0) return;
1960 char branchname[10];
1961 sprintf(branchname,"%s2",GetName());
1963 // Get the pointer to the header
1964 const char *cH = strstr(option,"H");
1966 if (fTrackHits && gAlice->TreeH() && cH) {
1967 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHits,
1968 gAlice->TreeH(),fBufferSize,1);
1969 gAlice->TreeH()->GetListOfBranches()->Add(branch);
1970 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
1972 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
1973 TDirectory *wd = gDirectory;
1975 TIter next( b->GetListOfBranches());
1976 while ((b=(TBranch*)next())) {
1980 cout << "Diverting branch " << branchname << " to file " << file << endl;
1985 void AliTPC::SetTreeAddress()
1987 if (fHitType&1) AliDetector::SetTreeAddress();
1988 if (fHitType&2) SetTreeAddress2();
1991 void AliTPC::SetTreeAddress2()
1994 // Set branch address for the TrackHits Tree
1997 char branchname[20];
1998 sprintf(branchname,"%s2",GetName());
2000 // Branch address for hit tree
2001 TTree *treeH = gAlice->TreeH();
2003 branch = treeH->GetBranch(branchname);
2004 if (branch) branch->SetAddress(&fTrackHits);
2008 void AliTPC::FinishPrimary()
2010 if (fTrackHits) fTrackHits->FlushHitStack();
2014 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2017 // add hit to the list
2020 int primary = gAlice->GetPrimary(track);
2021 gAlice->Particle(primary)->SetBit(kKeepBit);
2025 gAlice->FlagTrack(track);
2027 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2028 //if (hit->fTrack!=rtrack)
2029 // cout<<"bad track number\n";
2031 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2032 hits[1],hits[2],(Int_t)hits[3]);
2035 void AliTPC::ResetHits()
2037 if (fHitType&1) AliDetector::ResetHits();
2038 if (fHitType&2) ResetHits2();
2041 void AliTPC::ResetHits2()
2045 if (fTrackHits) fTrackHits->Clear();
2048 AliHit* AliTPC::FirstHit(Int_t track)
2050 if (fHitType&2) return FirstHit2(track);
2051 return AliDetector::FirstHit(track);
2053 AliHit* AliTPC::NextHit()
2055 if (fHitType&2) return NextHit2();
2056 return AliDetector::NextHit();
2059 AliHit* AliTPC::FirstHit2(Int_t track)
2062 // Initialise the hit iterator
2063 // Return the address of the first hit for track
2064 // If track>=0 the track is read from disk
2065 // while if track<0 the first hit of the current
2066 // track is returned
2069 gAlice->ResetHits();
2070 gAlice->TreeH()->GetEvent(track);
2074 fTrackHits->First();
2075 return fTrackHits->GetHit();
2080 AliHit* AliTPC::NextHit2()
2083 //Return the next hit for the current track
2087 return fTrackHits->GetHit();
2093 void AliTPC::LoadPoints(Int_t)
2097 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2100 if(fHitType==1) AliDetector::LoadPoints(a);
2101 else LoadPoints2(a);
2108 void AliTPC::RemapTrackHitIDs(Int_t *map)
2110 if (!fTrackHits) return;
2111 AliObjectArray * arr = fTrackHits->fTrackHitsInfo;
2112 for (UInt_t i=0;i<arr->GetSize();i++){
2113 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2114 info->fTrackID = map[info->fTrackID];
2120 //_____________________________________________________________________________
2121 void AliTPC::LoadPoints2(Int_t)
2124 // Store x, y, z of all hits in memory
2126 if (fTrackHits == 0) return;
2128 Int_t nhits = fTrackHits->GetEntriesFast();
2129 if (nhits == 0) return;
2130 Int_t tracks = gAlice->GetNtrack();
2131 if (fPoints == 0) fPoints = new TObjArray(tracks);
2134 Int_t *ntrk=new Int_t[tracks];
2135 Int_t *limi=new Int_t[tracks];
2136 Float_t **coor=new Float_t*[tracks];
2137 for(Int_t i=0;i<tracks;i++) {
2143 AliPoints *points = 0;
2146 Int_t chunk=nhits/4+1;
2148 // Loop over all the hits and store their position
2150 ahit = FirstHit2(-1);
2151 //for (Int_t hit=0;hit<nhits;hit++) {
2153 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2154 trk=ahit->GetTrack();
2155 if(ntrk[trk]==limi[trk]) {
2157 // Initialise a new track
2158 fp=new Float_t[3*(limi[trk]+chunk)];
2160 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2161 delete [] coor[trk];
2168 fp[3*ntrk[trk] ] = ahit->X();
2169 fp[3*ntrk[trk]+1] = ahit->Y();
2170 fp[3*ntrk[trk]+2] = ahit->Z();
2175 for(trk=0; trk<tracks; ++trk) {
2177 points = new AliPoints();
2178 points->SetMarkerColor(GetMarkerColor());
2179 points->SetMarkerSize(GetMarkerSize());
2180 points->SetDetector(this);
2181 points->SetParticle(trk);
2182 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2183 fPoints->AddAt(points,trk);
2184 delete [] coor[trk];
2194 //_____________________________________________________________________________
2195 void AliTPC::LoadPoints3(Int_t)
2198 // Store x, y, z of all hits in memory
2199 // - only intersection point with pad row
2200 if (fTrackHits == 0) return;
2202 Int_t nhits = fTrackHits->GetEntriesFast();
2203 if (nhits == 0) return;
2204 Int_t tracks = gAlice->GetNtrack();
2205 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2206 fPoints->Expand(2*tracks);
2209 Int_t *ntrk=new Int_t[tracks];
2210 Int_t *limi=new Int_t[tracks];
2211 Float_t **coor=new Float_t*[tracks];
2212 for(Int_t i=0;i<tracks;i++) {
2218 AliPoints *points = 0;
2221 Int_t chunk=nhits/4+1;
2223 // Loop over all the hits and store their position
2225 ahit = FirstHit2(-1);
2226 //for (Int_t hit=0;hit<nhits;hit++) {
2230 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2231 trk=ahit->GetTrack();
2232 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2233 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2234 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2235 if (currentrow!=lastrow){
2236 lastrow = currentrow;
2237 //later calculate intersection point
2238 if(ntrk[trk]==limi[trk]) {
2240 // Initialise a new track
2241 fp=new Float_t[3*(limi[trk]+chunk)];
2243 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2244 delete [] coor[trk];
2251 fp[3*ntrk[trk] ] = ahit->X();
2252 fp[3*ntrk[trk]+1] = ahit->Y();
2253 fp[3*ntrk[trk]+2] = ahit->Z();
2260 for(trk=0; trk<tracks; ++trk) {
2262 points = new AliPoints();
2263 points->SetMarkerColor(GetMarkerColor()+1);
2264 points->SetMarkerStyle(5);
2265 points->SetMarkerSize(0.2);
2266 points->SetDetector(this);
2267 points->SetParticle(trk);
2268 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2269 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2270 fPoints->AddAt(points,tracks+trk);
2271 delete [] coor[trk];
2282 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2286 //fill clones array with intersection of current point with the
2291 const Int_t kcmaxhits=30000;
2292 TVector * xxxx = new TVector(kcmaxhits*4);
2293 TVector & xxx = *xxxx;
2294 Int_t maxhits = kcmaxhits;
2297 AliTPChit * tpcHit=0;
2298 tpcHit = (AliTPChit*)FirstHit2(-1);
2299 Int_t currentIndex=0;
2300 Int_t lastrow=-1; //last writen row
2303 if (tpcHit==0) continue;
2304 sector=tpcHit->fSector; // sector number
2305 ipart=tpcHit->Track();
2309 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2310 Int_t index[3]={1,sector,0};
2311 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2312 if (currentrow<0) continue;
2313 if (lastrow<0) lastrow=currentrow;
2314 if (currentrow==lastrow){
2315 if ( currentIndex>=maxhits){
2317 xxx.ResizeTo(4*maxhits);
2319 xxx(currentIndex*4)=x[0];
2320 xxx(currentIndex*4+1)=x[1];
2321 xxx(currentIndex*4+2)=x[2];
2322 xxx(currentIndex*4+3)=tpcHit->fQ;
2326 if (currentIndex>2){
2338 for (Int_t index=0;index<currentIndex;index++){
2340 x=x2=x3=x4=xxx(index*4);
2348 sumy+=xxx(index*4+1);
2349 sumxy+=xxx(index*4+1)*x;
2350 sumx2y+=xxx(index*4+1)*x2;
2351 sumz+=xxx(index*4+2);
2352 sumxz+=xxx(index*4+2)*x;
2353 sumx2z+=xxx(index*4+2)*x2;
2354 sumq+=xxx(index*4+3);
2356 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2357 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2358 sumx2*(sumx*sumx3-sumx2*sumx2);
2360 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2361 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2362 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2363 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2365 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2366 sumx2*(sumx*sumx2y-sumx2*sumxy);
2367 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2368 sumx2*(sumx*sumx2z-sumx2*sumxz);
2370 Float_t y=detay/det+centralPad;
2371 Float_t z=detaz/det;
2372 Float_t by=detby/det; //y angle
2373 Float_t bz=detbz/det; //z angle
2374 sumy/=Float_t(currentIndex);
2375 sumz/=Float_t(currentIndex);
2377 AliComplexCluster cl;
2383 cl.fTracks[0]=ipart;
2385 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2386 if (row!=0) row->InsertCluster(&cl);
2389 } //end of calculating cluster for given row
2391 } // end of loop over hits