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.33 2001/04/03 12:40:43 kowal2
21 Revision 1.32 2001/03/12 17:47:36 hristov
22 Changes needed on Sun with CC 5.0
24 Revision 1.31 2001/03/12 08:21:50 kowal2
25 Corrected C++ bug in the material definitions
27 Revision 1.30 2001/03/01 17:34:47 kowal2
28 Correction due to the accuracy problem
30 Revision 1.29 2001/02/28 16:34:40 kowal2
31 Protection against nonphysical values of the avalanche size,
34 Revision 1.28 2001/01/26 19:57:19 hristov
35 Major upgrade of AliRoot code
37 Revision 1.27 2001/01/13 17:29:33 kowal2
38 Sun compiler correction
40 Revision 1.26 2001/01/10 07:59:43 kowal2
41 Corrections to load points from the noncompressed hits.
43 Revision 1.25 2000/11/02 07:25:31 kowal2
44 Changes due to the new hit structure.
47 Revision 1.24 2000/10/05 16:06:09 kowal2
48 Forward declarations. Changes due to a new class AliComplexCluster.
50 Revision 1.23 2000/10/02 21:28:18 fca
51 Removal of useless dependecies via forward declarations
53 Revision 1.22 2000/07/10 20:57:39 hristov
54 Update of TPC code and macros by M.Kowalski
56 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
57 Changes to obey the coding rules
59 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
60 Splitted from AliTPCtracking
62 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
63 Changed parameter settings
65 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
67 Defaults loaded automatically (hard-wired)
68 Optional parameters can be set via macro called in the constructor
70 Revision 1.19 2000/04/18 19:00:59 fca
71 Small bug fixes to TPC files
73 Revision 1.18 2000/04/17 09:37:33 kowal2
74 removed obsolete AliTPCDigitsDisplay.C
76 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
78 New, experimental data structure from M. Ivanov
79 New tracking algorithm
80 Different pad geometry for different sectors
81 Digitization rewritten
83 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
84 Not used anymore - removed
86 Revision 1.17 2000/01/19 17:17:30 fca
87 Introducing a list of lists of hits -- more hits allowed for detector now
89 Revision 1.16 1999/11/05 09:29:23 fca
90 Accept only signals > 0
92 Revision 1.15 1999/10/08 06:26:53 fca
93 Removed ClustersIndex - not used anymore
95 Revision 1.14 1999/09/29 09:24:33 fca
96 Introduction of the Copyright and cvs Log
100 ///////////////////////////////////////////////////////////////////////////////
102 // Time Projection Chamber //
103 // This class contains the basic functions for the Time Projection Chamber //
104 // detector. Functions specific to one particular geometry are //
105 // contained in the derived classes //
109 <img src="picts/AliTPCClass.gif">
114 ///////////////////////////////////////////////////////////////////////////////
122 #include <TGeometry.h>
125 #include <TObjectTable.h>
126 #include "TParticle.h"
130 #include <iostream.h>
137 #include "AliTPCParamSR.h"
138 #include "AliTPCPRF2D.h"
139 #include "AliTPCRF1D.h"
140 #include "AliDigits.h"
141 #include "AliSimDigits.h"
142 #include "AliTPCTrackHits.h"
143 #include "AliPoints.h"
144 #include "AliArrayBranch.h"
147 #include "AliTPCDigitsArray.h"
148 #include "AliComplexCluster.h"
149 #include "AliClusters.h"
150 #include "AliTPCClustersRow.h"
151 #include "AliTPCClustersArray.h"
153 #include "AliTPCcluster.h"
154 #include "AliTPCclusterer.h"
155 #include "AliTPCtracker.h"
157 #include <TInterpreter.h>
164 //_____________________________________________________________________________
168 // Default constructor
183 //_____________________________________________________________________________
184 AliTPC::AliTPC(const char *name, const char *title)
185 : AliDetector(name,title)
188 // Standard constructor
192 // Initialise arrays of hits and digits
193 fHits = new TClonesArray("AliTPChit", 176);
194 gAlice->AddHitList(fHits);
199 fTrackHits = new AliTPCTrackHits; //MI - 13.09.2000
200 fTrackHits->SetHitPrecision(0.002);
201 fTrackHits->SetStepPrecision(0.003);
202 fTrackHits->SetMaxDistance(100);
206 // Initialise counters
212 // Initialise color attributes
213 SetMarkerColor(kYellow);
216 // Set TPC parameters
219 if (!strcmp(title,"Default")) {
220 fTPCParam = new AliTPCParamSR;
222 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
228 //_____________________________________________________________________________
239 delete fTrackHits; //MI 15.09.2000
242 //_____________________________________________________________________________
243 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
246 // Add a hit to the list
248 // TClonesArray &lhits = *fHits;
249 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
251 TClonesArray &lhits = *fHits;
252 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
255 AddHit2(track,vol,hits);
258 //_____________________________________________________________________________
259 void AliTPC::BuildGeometry()
263 // Build TPC ROOT TNode geometry for the event display
268 const int kColorTPC=19;
269 char name[5], title[25];
270 const Double_t kDegrad=TMath::Pi()/180;
271 const Double_t kRaddeg=180./TMath::Pi();
274 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
275 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
277 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
278 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
280 Int_t nLo = fTPCParam->GetNInnerSector()/2;
281 Int_t nHi = fTPCParam->GetNOuterSector()/2;
283 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
284 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
285 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
286 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
289 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
290 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
296 // Get ALICE top node
299 nTop=gAlice->GetGeometry()->GetNode("alice");
303 rl = fTPCParam->GetInnerRadiusLow();
304 ru = fTPCParam->GetInnerRadiusUp();
308 sprintf(name,"LS%2.2d",i);
310 sprintf(title,"TPC low sector %3d",i);
313 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
314 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
315 tubs->SetNumberOfDivisions(1);
317 nNode = new TNode(name,title,name,0,0,0,"");
318 nNode->SetLineColor(kColorTPC);
324 rl = fTPCParam->GetOuterRadiusLow();
325 ru = fTPCParam->GetOuterRadiusUp();
328 sprintf(name,"US%2.2d",i);
330 sprintf(title,"TPC upper sector %d",i);
332 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
333 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
334 tubs->SetNumberOfDivisions(1);
336 nNode = new TNode(name,title,name,0,0,0,"");
337 nNode->SetLineColor(kColorTPC);
343 //_____________________________________________________________________________
344 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
347 // Calculate distance from TPC to mouse on the display
353 void AliTPC::Clusters2Tracks(TFile *of) {
354 //-----------------------------------------------------------------
355 // This is a track finder.
356 //-----------------------------------------------------------------
357 AliTPCtracker tracker(fTPCParam);
358 tracker.Clusters2Tracks(gFile,of);
361 //_____________________________________________________________________________
362 void AliTPC::CreateMaterials()
364 //-----------------------------------------------
365 // Create Materials for for TPC simulations
366 //-----------------------------------------------
368 //-----------------------------------------------------------------
369 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
370 //-----------------------------------------------------------------
372 Int_t iSXFLD=gAlice->Field()->Integ();
373 Float_t sXMGMX=gAlice->Field()->Max();
375 Float_t amat[5]; // atomic numbers
376 Float_t zmat[5]; // z
377 Float_t wmat[5]; // proportions
383 //***************** Gases *************************
385 //-------------------------------------------------
387 //-------------------------------------------------
398 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
408 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
411 //--------------------------------------------------------------
413 //--------------------------------------------------------------
430 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
432 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
447 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
449 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
465 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
467 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
469 //----------------------------------------------------------------
470 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
471 //----------------------------------------------------------------
477 Float_t rho,absl,X0,buf[1];
481 for(nc = 0;nc<fNoComp;nc++)
484 // retrive material constants
486 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
491 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
493 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
494 density += fMixtProp[nc]*rho; // density of the mixture
498 // mixture proportions by weight!
500 for(nc = 0;nc<fNoComp;nc++)
503 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
505 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
506 apure[nnc] : amol[nnc])/am;
510 // Drift gases 1 - nonsensitive, 2 - sensitive
512 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
513 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
522 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
525 //----------------------------------------------------------------------
527 //----------------------------------------------------------------------
549 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
571 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
589 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
607 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
609 // G10 60% SiO2 + 40% epoxy, I use A and Z for SiO2
622 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz
624 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
626 AliMaterial(39,"G10",amat[0],zmat[0],density,999.,999.);
635 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
644 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
653 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
671 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
690 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
692 // Lexan - I use it for epoxy - similar X0 C16H14O3
709 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
711 //----------------------------------------------------------
712 // tracking media for gases
713 //----------------------------------------------------------
715 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
716 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
717 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
718 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
720 //-----------------------------------------------------------
721 // tracking media for solids
722 //-----------------------------------------------------------
724 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
725 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
726 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
727 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
728 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
729 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
730 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
731 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
732 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
733 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
734 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
739 void AliTPC::Digits2Clusters(TFile *of)
741 //-----------------------------------------------------------------
742 // This is a simple cluster finder.
743 //-----------------------------------------------------------------
744 AliTPCclusterer::Digits2Clusters(fTPCParam,of);
747 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
748 extern Double_t SigmaZ2(Double_t, Double_t);
749 //_____________________________________________________________________________
750 void AliTPC::Hits2Clusters(TFile *of)
752 //--------------------------------------------------------
753 // TPC simple cluster generator from hits
754 // obtained from the TPC Fast Simulator
755 // The point errors are taken from the parametrization
756 //--------------------------------------------------------
758 //-----------------------------------------------------------------
759 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
760 //-----------------------------------------------------------------
761 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
762 // Jouri.Belikov@cern.ch
763 //----------------------------------------------------------------
765 /////////////////////////////////////////////////////////////////////////////
767 //---------------------------------------------------------------------
768 // ALICE TPC Cluster Parameters
769 //--------------------------------------------------------------------
773 // Cluster width in rphi
774 const Float_t kACrphi=0.18322;
775 const Float_t kBCrphi=0.59551e-3;
776 const Float_t kCCrphi=0.60952e-1;
777 // Cluster width in z
778 const Float_t kACz=0.19081;
779 const Float_t kBCz=0.55938e-3;
780 const Float_t kCCz=0.30428;
782 TDirectory *savedir=gDirectory;
785 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
790 printf("AliTPCParam MUST be created firstly\n");
794 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
796 TParticle *particle; // pointer to a given particle
797 AliTPChit *tpcHit; // pointer to a sigle TPC hit
801 Float_t pl,pt,tanth,rpad,ratio;
804 //---------------------------------------------------------------
805 // Get the access to the tracks
806 //---------------------------------------------------------------
808 TTree *tH = gAlice->TreeH();
809 Stat_t ntracks = tH->GetEntries();
811 //Switch to the output file
814 fTPCParam->Write(fTPCParam->GetTitle());
815 AliTPCClustersArray carray;
816 carray.Setup(fTPCParam);
817 carray.SetClusterType("AliTPCcluster");
820 Int_t nclusters=0; //cluster counter
822 //------------------------------------------------------------
823 // Loop over all sectors (72 sectors for 20 deg
824 // segmentation for both lower and upper sectors)
825 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
826 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
828 // First cluster for sector 0 starts at "0"
829 //------------------------------------------------------------
831 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
833 fTPCParam->AdjustCosSin(isec,cph,sph);
835 //------------------------------------------------------------
837 //------------------------------------------------------------
839 for(Int_t track=0;track<ntracks;track++){
843 // Get number of the TPC hits
845 // nhits=fHits->GetEntriesFast();
848 tpcHit = (AliTPChit*)FirstHit(-1);
852 // for(Int_t hit=0;hit<nhits;hit++){
853 //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
857 if (tpcHit->fQ == 0.) {
858 tpcHit = (AliTPChit*) NextHit();
859 continue; //information about track (I.Belikov)
861 sector=tpcHit->fSector; // sector number
864 // if(sector != isec) continue; //terminate iteration
867 tpcHit = (AliTPChit*) NextHit();
870 ipart=tpcHit->Track();
871 particle=gAlice->Particle(ipart);
874 if(pt < 1.e-9) pt=1.e-9;
876 tanth = TMath::Abs(tanth);
877 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
878 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
880 // space-point resolutions
882 sigmaRphi=SigmaY2(rpad,tanth,pt);
883 sigmaZ =SigmaZ2(rpad,tanth );
887 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
888 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
890 // temporary protection
892 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
893 if(sigmaZ < 0.) sigmaZ=0.4e-3;
894 if(clRphi < 0.) clRphi=2.5e-3;
895 if(clZ < 0.) clZ=2.5e-5;
900 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
901 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
903 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
904 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
905 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
906 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
907 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
908 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
909 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
910 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
911 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
912 xyz[2]=sigmaRphi; // fSigmaY2
913 xyz[3]=sigmaZ; // fSigmaZ2
914 xyz[4]=tpcHit->fQ; // q
916 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
917 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
919 Int_t tracks[3]={tpcHit->Track(), -1, -1};
920 AliTPCcluster cluster(tracks,xyz);
922 clrow->InsertCluster(&cluster); nclusters++;
924 tpcHit = (AliTPChit*)NextHit();
927 } // end of loop over hits
929 } // end of loop over tracks
931 Int_t nrows=fTPCParam->GetNRow(isec);
932 for (Int_t irow=0; irow<nrows; irow++) {
933 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
934 if (!clrow) continue;
935 carray.StoreRow(isec,irow);
936 carray.ClearRow(isec,irow);
939 } // end of loop over sectors
941 cerr<<"Number of made clusters : "<<nclusters<<" \n";
943 carray.GetTree()->Write();
945 savedir->cd(); //switch back to the input file
949 //_________________________________________________________________
950 void AliTPC::Hits2ExactClustersSector(Int_t isec)
952 //--------------------------------------------------------
953 //calculate exact cross point of track and given pad row
954 //resulting values are expressed in "digit" coordinata
955 //--------------------------------------------------------
957 //-----------------------------------------------------------------
958 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
959 //-----------------------------------------------------------------
961 if (fClustersArray==0){
965 TParticle *particle; // pointer to a given particle
966 AliTPChit *tpcHit; // pointer to a sigle TPC hit
969 const Int_t kcmaxhits=30000;
970 TVector * xxxx = new TVector(kcmaxhits*4);
971 TVector & xxx = *xxxx;
972 Int_t maxhits = kcmaxhits;
973 //construct array for each padrow
974 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
975 fClustersArray->CreateRow(isec,i);
977 //---------------------------------------------------------------
978 // Get the access to the tracks
979 //---------------------------------------------------------------
981 TTree *tH = gAlice->TreeH();
982 Stat_t ntracks = tH->GetEntries();
983 Int_t npart = gAlice->GetNtrack();
985 //------------------------------------------------------------
987 //------------------------------------------------------------
989 for(Int_t track=0;track<ntracks;track++){
993 // Get number of the TPC hits and a pointer
996 nhits=fHits->GetEntriesFast();
1000 Int_t currentIndex=0;
1001 Int_t lastrow=-1; //last writen row
1002 for(Int_t hit=0;hit<nhits;hit++){
1003 tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
1004 if (tpcHit==0) continue;
1005 sector=tpcHit->fSector; // sector number
1006 if(sector != isec) continue;
1007 ipart=tpcHit->Track();
1008 if (ipart<npart) particle=gAlice->Particle(ipart);
1012 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1013 Int_t index[3]={1,isec,0};
1014 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1015 if (currentrow<0) continue;
1016 if (lastrow<0) lastrow=currentrow;
1017 if (currentrow==lastrow){
1018 if ( currentIndex>=maxhits){
1020 xxx.ResizeTo(4*maxhits);
1022 xxx(currentIndex*4)=x[0];
1023 xxx(currentIndex*4+1)=x[1];
1024 xxx(currentIndex*4+2)=x[2];
1025 xxx(currentIndex*4+3)=tpcHit->fQ;
1029 if (currentIndex>2){
1041 for (Int_t index=0;index<currentIndex;index++){
1043 x=x2=x3=x4=xxx(index*4);
1051 sumy+=xxx(index*4+1);
1052 sumxy+=xxx(index*4+1)*x;
1053 sumx2y+=xxx(index*4+1)*x2;
1054 sumz+=xxx(index*4+2);
1055 sumxz+=xxx(index*4+2)*x;
1056 sumx2z+=xxx(index*4+2)*x2;
1057 sumq+=xxx(index*4+3);
1059 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1060 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1061 sumx2*(sumx*sumx3-sumx2*sumx2);
1063 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1064 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1065 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1066 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1068 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1069 sumx2*(sumx*sumx2y-sumx2*sumxy);
1070 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1071 sumx2*(sumx*sumx2z-sumx2*sumxz);
1073 Float_t y=detay/det+centralPad;
1074 Float_t z=detaz/det;
1075 Float_t by=detby/det; //y angle
1076 Float_t bz=detbz/det; //z angle
1077 sumy/=Float_t(currentIndex);
1078 sumz/=Float_t(currentIndex);
1079 AliComplexCluster cl;
1085 cl.fTracks[0]=ipart;
1087 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1088 if (row!=0) row->InsertCluster(&cl);
1091 } //end of calculating cluster for given row
1095 } // end of loop over hits
1096 } // end of loop over tracks
1097 //write padrows to tree
1098 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1099 fClustersArray->StoreRow(isec,ii);
1100 fClustersArray->ClearRow(isec,ii);
1105 //___________________________________________
1106 void AliTPC::SDigits2Digits()
1108 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1109 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1110 AliTPCPRF2D * prfouter = new AliTPCPRF2D;
1111 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1113 TDirectory *cwd = gDirectory;
1114 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1115 rf->SetOffset(3*param->GetZSigma());
1117 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1119 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n";
1122 prfinner->Read("prf_07504_Gati_056068_d02");
1123 prfouter->Read("prf_10006_Gati_047051_d03");
1127 param->SetInnerPRF(prfinner);
1128 param->SetOuterPRF(prfouter);
1129 param->SetTimeRF(rf);
1133 cerr<<"Digitizing TPC...\n";
1135 //setup TPCDigitsArray
1136 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1137 arr->SetClass("AliSimDigits");
1141 arr->MakeTree(fDigitsFile);
1143 SetDigitsArray(arr);
1147 // Hits2DigitsSector(1);
1148 // Hits2DigitsSector(2);
1149 // Hits2DigitsSector(3);
1150 // Hits2DigitsSector(1+18);
1151 // Hits2DigitsSector(2+18);
1152 // Hits2DigitsSector(3+18);
1154 // Hits2DigitsSector(36+1);
1155 // Hits2DigitsSector(36+2);
1156 // Hits2DigitsSector(36+3);
1157 // Hits2DigitsSector(36+1+18);
1158 // Hits2DigitsSector(36+2+18);
1159 // Hits2DigitsSector(36+3+18);
1164 sprintf(treeName,"TreeD_%s",param->GetTitle());
1165 GetDigitsArray()->GetTree()->Write(treeName,TObject::kOverwrite);
1168 //__________________________________________________________________
1169 void AliTPC::Hits2Digits()
1171 //----------------------------------------------------
1172 // Loop over all sectors
1173 //----------------------------------------------------
1176 printf("AliTPCParam MUST be created firstly\n");
1180 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1185 //_____________________________________________________________________________
1186 void AliTPC::Hits2DigitsSector(Int_t isec)
1188 //-------------------------------------------------------------------
1189 // TPC conversion from hits to digits.
1190 //-------------------------------------------------------------------
1192 //-----------------------------------------------------------------
1193 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1194 //-----------------------------------------------------------------
1196 //-------------------------------------------------------
1197 // Get the access to the track hits
1198 //-------------------------------------------------------
1201 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1202 Stat_t ntracks = tH->GetEntries();
1206 //-------------------------------------------
1207 // Only if there are any tracks...
1208 //-------------------------------------------
1212 //printf("*** Processing sector number %d ***\n",isec);
1214 Int_t nrows =fTPCParam->GetNRow(isec);
1216 row= new TObjArray* [nrows];
1218 MakeSector(isec,nrows,tH,ntracks,row);
1220 //--------------------------------------------------------
1221 // Digitize this sector, row by row
1222 // row[i] is the pointer to the TObjArray of TVectors,
1223 // each one containing electrons accepted on this
1224 // row, assigned into tracks
1225 //--------------------------------------------------------
1229 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1231 for (i=0;i<nrows;i++){
1233 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1235 DigitizeRow(i,isec,row);
1237 fDigitsArray->StoreRow(isec,i);
1239 Int_t ndig = dig->GetDigitSize();
1242 //printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1244 fDigitsArray->ClearRow(isec,i);
1247 } // end of the sector digitization
1249 for(i=0;i<nrows;i++){
1254 delete [] row; // delete the array of pointers to TObjArray-s
1258 } // end of Hits2DigitsSector
1261 //_____________________________________________________________________________
1262 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1264 //-----------------------------------------------------------
1265 // Single row digitization, coupling from the neighbouring
1266 // rows taken into account
1267 //-----------------------------------------------------------
1269 //-----------------------------------------------------------------
1270 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1271 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1272 //-----------------------------------------------------------------
1275 Float_t zerosup = fTPCParam->GetZeroSup();
1276 Int_t nrows =fTPCParam->GetNRow(isec);
1277 fCurrentIndex[1]= isec;
1280 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1281 Int_t nofTbins = fTPCParam->GetMaxTBin();
1282 Int_t indexRange[4];
1284 // Integrated signal for this row
1285 // and a single track signal
1287 TMatrix *m1 = new TMatrix(0,nofPads,0,nofTbins); // integrated
1288 TMatrix *m2 = new TMatrix(0,nofPads,0,nofTbins); // single
1290 TMatrix &total = *m1;
1292 // Array of pointers to the label-signal list
1294 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1295 Float_t **pList = new Float_t* [nofDigits];
1299 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1303 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1304 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1305 for (Int_t row= row1;row<=row2;row++){
1306 Int_t nTracks= rows[row]->GetEntries();
1307 for (i1=0;i1<nTracks;i1++){
1308 fCurrentIndex[2]= row;
1309 fCurrentIndex[3]=irow;
1311 m2->Zero(); // clear single track signal matrix
1312 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1313 GetList(trackLabel,nofPads,m2,indexRange,pList);
1315 else GetSignal(rows[row],i1,0,m1,indexRange);
1321 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1322 for(Int_t ip=0;ip<nofPads;ip++){
1323 for(Int_t it=0;it<nofTbins;it++){
1325 Float_t q = total(ip,it);
1327 Int_t gi =it*nofPads+ip; // global index
1329 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
1333 if(q <=zerosup) continue; // do not fill zeros
1334 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1337 // "real" signal or electronic noise (list = -1)?
1340 for(Int_t j1=0;j1<3;j1++){
1341 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1346 <A NAME="AliDigits"></A>
1347 using of AliDigits object
1350 dig->SetDigitFast((Short_t)q,it,ip);
1351 if (fDigitsArray->IsSimulated())
1353 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1354 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1355 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1359 } // end of loop over time buckets
1360 } // end of lop over pads
1363 // This row has been digitized, delete nonused stuff
1366 for(lp=0;lp<nofDigits;lp++){
1367 if(pList[lp]) delete [] pList[lp];
1376 } // end of DigitizeRow
1378 //_____________________________________________________________________________
1380 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, TMatrix *m1, TMatrix *m2,
1384 //---------------------------------------------------------------
1385 // Calculates 2-D signal (pad,time) for a single track,
1386 // returns a pointer to the signal matrix and the track label
1387 // No digitization is performed at this level!!!
1388 //---------------------------------------------------------------
1390 //-----------------------------------------------------------------
1391 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1392 // Modified: Marian Ivanov
1393 //-----------------------------------------------------------------
1397 tv = (TVector*)p1->At(ntr); // pointer to a track
1400 Float_t label = v(0);
1401 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1403 Int_t nElectrons = (tv->GetNrows()-1)/4;
1404 indexRange[0]=9999; // min pad
1405 indexRange[1]=-1; // max pad
1406 indexRange[2]=9999; //min time
1407 indexRange[3]=-1; // max time
1409 // Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1411 TMatrix &signal = *m1;
1412 TMatrix &total = *m2;
1414 // Loop over all electrons
1416 for(Int_t nel=0; nel<nElectrons; nel++){
1418 Float_t aval = v(idx+4);
1419 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1420 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1421 Int_t n = fTPCParam->CalcResponse(xyz,fCurrentIndex,fCurrentIndex[3]);
1423 if (n>0) for (Int_t i =0; i<n; i++){
1424 Int_t *index = fTPCParam->GetResBin(i);
1425 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1426 if ( ( pad<(fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]))) && (pad>0)) {
1427 Int_t time=index[2];
1428 Float_t weight = fTPCParam->GetResWeight(i); //we normalise response to ADC channel
1429 weight *= eltoadcfac;
1431 if (m1!=0) signal(pad,time)+=weight;
1432 total(pad,time)+=weight;
1433 indexRange[0]=TMath::Min(indexRange[0],pad);
1434 indexRange[1]=TMath::Max(indexRange[1],pad);
1435 indexRange[2]=TMath::Min(indexRange[2],time);
1436 indexRange[3]=TMath::Max(indexRange[3],time);
1439 } // end of loop over electrons
1441 return label; // returns track label when finished
1444 //_____________________________________________________________________________
1445 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *indexRange,
1448 //----------------------------------------------------------------------
1449 // Updates the list of tracks contributing to digits for a given row
1450 //----------------------------------------------------------------------
1452 //-----------------------------------------------------------------
1453 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1454 //-----------------------------------------------------------------
1456 TMatrix &signal = *m;
1458 // lop over nonzero digits
1460 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1461 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1464 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1466 if(signal(ip,it)<0.5) continue;
1469 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1471 if(!pList[globalIndex]){
1474 // Create new list (6 elements - 3 signals and 3 labels),
1477 pList[globalIndex] = new Float_t [6];
1481 *pList[globalIndex] = -1.;
1482 *(pList[globalIndex]+1) = -1.;
1483 *(pList[globalIndex]+2) = -1.;
1484 *(pList[globalIndex]+3) = -1.;
1485 *(pList[globalIndex]+4) = -1.;
1486 *(pList[globalIndex]+5) = -1.;
1489 *pList[globalIndex] = label;
1490 *(pList[globalIndex]+3) = signal(ip,it);
1494 // check the signal magnitude
1496 Float_t highest = *(pList[globalIndex]+3);
1497 Float_t middle = *(pList[globalIndex]+4);
1498 Float_t lowest = *(pList[globalIndex]+5);
1501 // compare the new signal with already existing list
1504 if(signal(ip,it)<lowest) continue; // neglect this track
1508 if (signal(ip,it)>highest){
1509 *(pList[globalIndex]+5) = middle;
1510 *(pList[globalIndex]+4) = highest;
1511 *(pList[globalIndex]+3) = signal(ip,it);
1513 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1514 *(pList[globalIndex]+1) = *pList[globalIndex];
1515 *pList[globalIndex] = label;
1517 else if (signal(ip,it)>middle){
1518 *(pList[globalIndex]+5) = middle;
1519 *(pList[globalIndex]+4) = signal(ip,it);
1521 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1522 *(pList[globalIndex]+1) = label;
1525 *(pList[globalIndex]+5) = signal(ip,it);
1526 *(pList[globalIndex]+2) = label;
1530 } // end of loop over pads
1531 } // end of loop over time bins
1536 //___________________________________________________________________
1537 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1538 Stat_t ntracks,TObjArray **row)
1541 //-----------------------------------------------------------------
1542 // Prepares the sector digitization, creates the vectors of
1543 // tracks for each row of this sector. The track vector
1544 // contains the track label and the position of electrons.
1545 //-----------------------------------------------------------------
1547 //-----------------------------------------------------------------
1548 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1549 //-----------------------------------------------------------------
1551 Float_t gasgain = fTPCParam->GetGasGain();
1555 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1558 if (fHitType&2) branch = TH->GetBranch("TPC2");
1559 else branch = TH->GetBranch("TPC");
1562 //----------------------------------------------
1563 // Create TObjArray-s, one for each row,
1564 // each TObjArray will store the TVectors
1565 // of electrons, one TVector per each track.
1566 //----------------------------------------------
1568 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
1569 TVector **tracks = new TVector* [nrows]; //pointers to the track vectors
1570 for(i=0; i<nrows; i++){
1571 row[i] = new TObjArray;
1578 //--------------------------------------------------------------------
1579 // Loop over tracks, the "track" contains the full history
1580 //--------------------------------------------------------------------
1582 Int_t previousTrack,currentTrack;
1583 previousTrack = -1; // nothing to store so far!
1585 for(Int_t track=0;track<ntracks;track++){
1586 Bool_t isInSector=kTRUE;
1591 TBranch * br = TH->GetBranch("fTrackHitsInfo");
1592 br->GetEvent(track);
1593 AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
1594 for (UInt_t j=0;j<ar->GetSize();j++){
1595 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
1598 if (!isInSector) continue;
1600 branch->GetEntry(track); // get next track
1604 tpcHit = (AliTPChit*)FirstHit(-1);
1606 //--------------------------------------------------------------
1608 //--------------------------------------------------------------
1613 Int_t sector=tpcHit->fSector; // sector number
1614 // if(sector != isec) continue;
1616 tpcHit = (AliTPChit*) NextHit();
1620 currentTrack = tpcHit->Track(); // track number
1623 if(currentTrack != previousTrack){
1625 // store already filled fTrack
1627 for(i=0;i<nrows;i++){
1628 if(previousTrack != -1){
1629 if(nofElectrons[i]>0){
1630 TVector &v = *tracks[i];
1631 v(0) = previousTrack;
1632 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1633 row[i]->Add(tracks[i]);
1636 delete tracks[i]; // delete empty TVector
1642 tracks[i] = new TVector(481); // TVectors for the next fTrack
1644 } // end of loop over rows
1646 previousTrack=currentTrack; // update track label
1649 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1651 //---------------------------------------------------
1652 // Calculate the electron attachment probability
1653 //---------------------------------------------------
1656 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
1657 /fTPCParam->GetDriftV();
1659 Float_t attProb = fTPCParam->GetAttCoef()*
1660 fTPCParam->GetOxyCont()*time; // fraction!
1662 //-----------------------------------------------
1663 // Loop over electrons
1664 //-----------------------------------------------
1667 for(Int_t nel=0;nel<qI;nel++){
1668 // skip if electron lost due to the attachment
1669 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1674 // protection for the nonphysical avalanche size (10**6 maximum)
1676 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
1677 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
1680 TransportElectron(xyz,index); //MI change -august
1682 fTPCParam->GetPadRow(xyz,index); //MI change august
1683 rowNumber = index[2];
1684 //transform position to local digit coordinates
1685 //relative to nearest pad row
1686 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
1687 nofElectrons[rowNumber]++;
1688 //----------------------------------
1689 // Expand vector if necessary
1690 //----------------------------------
1691 if(nofElectrons[rowNumber]>120){
1692 Int_t range = tracks[rowNumber]->GetNrows();
1693 if((nofElectrons[rowNumber])>(range-1)/4){
1695 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
1699 TVector &v = *tracks[rowNumber];
1700 Int_t idx = 4*nofElectrons[rowNumber]-3;
1702 v(idx)= xyz[0]; // X - pad row coordinate
1703 v(idx+1)=xyz[1]; // Y - pad coordinate (along the pad-row)
1704 v(idx+2)=xyz[2]; // Z - time bin coordinate
1705 v(idx+3)=xyz[3]; // avalanche size
1706 } // end of loop over electrons
1708 tpcHit = (AliTPChit*)NextHit();
1710 } // end of loop over hits
1711 } // end of loop over tracks
1714 // store remaining track (the last one) if not empty
1717 for(i=0;i<nrows;i++){
1718 if(nofElectrons[i]>0){
1719 TVector &v = *tracks[i];
1720 v(0) = previousTrack;
1721 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1722 row[i]->Add(tracks[i]);
1731 delete [] nofElectrons;
1734 } // end of MakeSector
1737 //_____________________________________________________________________________
1741 // Initialise TPC detector after definition of geometry
1746 for(i=0;i<35;i++) printf("*");
1747 printf(" TPC_INIT ");
1748 for(i=0;i<35;i++) printf("*");
1751 for(i=0;i<80;i++) printf("*");
1755 //_____________________________________________________________________________
1756 void AliTPC::MakeBranch(Option_t* option, char *file)
1759 // Create Tree branches for the TPC.
1761 Int_t buffersize = 4000;
1762 char branchname[10];
1763 sprintf(branchname,"%s",GetName());
1765 AliDetector::MakeBranch(option,file);
1767 const char *d = strstr(option,"D");
1769 if (fDigits && gAlice->TreeD() && d) {
1770 gAlice->MakeBranchInTree(gAlice->TreeD(),
1771 branchname, &fDigits, buffersize, file) ;
1774 if (fHitType&2) MakeBranch2(option,file); // MI change 14.09.2000
1777 //_____________________________________________________________________________
1778 void AliTPC::ResetDigits()
1781 // Reset number of digits and the digits array for this detector
1784 if (fDigits) fDigits->Clear();
1787 //_____________________________________________________________________________
1788 void AliTPC::SetSecAL(Int_t sec)
1790 //---------------------------------------------------
1791 // Activate/deactivate selection for lower sectors
1792 //---------------------------------------------------
1794 //-----------------------------------------------------------------
1795 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1796 //-----------------------------------------------------------------
1801 //_____________________________________________________________________________
1802 void AliTPC::SetSecAU(Int_t sec)
1804 //----------------------------------------------------
1805 // Activate/deactivate selection for upper sectors
1806 //---------------------------------------------------
1808 //-----------------------------------------------------------------
1809 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1810 //-----------------------------------------------------------------
1815 //_____________________________________________________________________________
1816 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
1818 //----------------------------------------
1819 // Select active lower sectors
1820 //----------------------------------------
1822 //-----------------------------------------------------------------
1823 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1824 //-----------------------------------------------------------------
1834 //_____________________________________________________________________________
1835 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
1836 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
1837 Int_t s11 , Int_t s12)
1839 //--------------------------------
1840 // Select active upper sectors
1841 //--------------------------------
1843 //-----------------------------------------------------------------
1844 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1845 //-----------------------------------------------------------------
1861 //_____________________________________________________________________________
1862 void AliTPC::SetSens(Int_t sens)
1865 //-------------------------------------------------------------
1866 // Activates/deactivates the sensitive strips at the center of
1867 // the pad row -- this is for the space-point resolution calculations
1868 //-------------------------------------------------------------
1870 //-----------------------------------------------------------------
1871 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1872 //-----------------------------------------------------------------
1878 void AliTPC::SetSide(Float_t side=0.)
1880 // choice of the TPC side
1885 //____________________________________________________________________________
1886 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
1887 Float_t p2,Float_t p3)
1890 // gax mixture definition
1904 //_____________________________________________________________________________
1906 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
1909 // electron transport taking into account:
1911 // 2.ExB at the wires
1912 // 3. nonisochronity
1914 // xyz and index must be already transformed to system 1
1917 fTPCParam->Transform1to2(xyz,index);
1920 Float_t driftl=xyz[2];
1921 if(driftl<0.01) driftl=0.01;
1922 driftl=TMath::Sqrt(driftl);
1923 Float_t sigT = driftl*(fTPCParam->GetDiffT());
1924 Float_t sigL = driftl*(fTPCParam->GetDiffL());
1925 xyz[0]=gRandom->Gaus(xyz[0],sigT);
1926 xyz[1]=gRandom->Gaus(xyz[1],sigT);
1927 xyz[2]=gRandom->Gaus(xyz[2],sigL);
1931 if (fTPCParam->GetMWPCReadout()==kTRUE){
1933 fTPCParam->Transform2to2NearestWire(xyz,index);
1934 Float_t dx=xyz[0]-x1;
1935 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
1937 //add nonisochronity (not implemented yet)
1941 ClassImp(AliTPCdigit)
1943 //_____________________________________________________________________________
1944 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
1948 // Creates a TPC digit object
1950 fSector = digits[0];
1951 fPadRow = digits[1];
1954 fSignal = digits[4];
1960 //_____________________________________________________________________________
1961 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
1965 // Creates a TPC hit object
1976 //________________________________________________________________________
1977 // Additional code because of the AliTPCTrackHits
1979 void AliTPC::MakeBranch2(Option_t *option,char *file)
1982 // Create a new branch in the current Root Tree
1983 // The branch of fHits is automatically split
1984 // MI change 14.09.2000
1985 if (fHitType&2==0) return;
1986 char branchname[10];
1987 sprintf(branchname,"%s2",GetName());
1989 // Get the pointer to the header
1990 const char *cH = strstr(option,"H");
1992 if (fTrackHits && gAlice->TreeH() && cH) {
1993 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHits,
1994 gAlice->TreeH(),fBufferSize,1);
1995 gAlice->TreeH()->GetListOfBranches()->Add(branch);
1996 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
1998 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
1999 TDirectory *wd = gDirectory;
2001 TIter next( b->GetListOfBranches());
2002 while ((b=(TBranch*)next())) {
2006 cout << "Diverting branch " << branchname << " to file " << file << endl;
2011 void AliTPC::SetTreeAddress()
2013 if (fHitType&1) AliDetector::SetTreeAddress();
2014 if (fHitType&2) SetTreeAddress2();
2017 void AliTPC::SetTreeAddress2()
2020 // Set branch address for the TrackHits Tree
2023 char branchname[20];
2024 sprintf(branchname,"%s2",GetName());
2026 // Branch address for hit tree
2027 TTree *treeH = gAlice->TreeH();
2029 branch = treeH->GetBranch(branchname);
2030 if (branch) branch->SetAddress(&fTrackHits);
2034 void AliTPC::FinishPrimary()
2036 if (fTrackHits) fTrackHits->FlushHitStack();
2040 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2043 // add hit to the list
2046 int primary = gAlice->GetPrimary(track);
2047 gAlice->Particle(primary)->SetBit(kKeepBit);
2051 gAlice->FlagTrack(track);
2053 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2054 //if (hit->fTrack!=rtrack)
2055 // cout<<"bad track number\n";
2057 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2058 hits[1],hits[2],(Int_t)hits[3]);
2061 void AliTPC::ResetHits()
2063 if (fHitType&1) AliDetector::ResetHits();
2064 if (fHitType&2) ResetHits2();
2067 void AliTPC::ResetHits2()
2071 if (fTrackHits) fTrackHits->Clear();
2074 AliHit* AliTPC::FirstHit(Int_t track)
2076 if (fHitType&2) return FirstHit2(track);
2077 return AliDetector::FirstHit(track);
2079 AliHit* AliTPC::NextHit()
2081 if (fHitType&2) return NextHit2();
2082 return AliDetector::NextHit();
2085 AliHit* AliTPC::FirstHit2(Int_t track)
2088 // Initialise the hit iterator
2089 // Return the address of the first hit for track
2090 // If track>=0 the track is read from disk
2091 // while if track<0 the first hit of the current
2092 // track is returned
2095 gAlice->ResetHits();
2096 gAlice->TreeH()->GetEvent(track);
2100 fTrackHits->First();
2101 return fTrackHits->GetHit();
2106 AliHit* AliTPC::NextHit2()
2109 //Return the next hit for the current track
2113 return fTrackHits->GetHit();
2119 void AliTPC::LoadPoints(Int_t)
2123 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2126 if(fHitType==1) AliDetector::LoadPoints(a);
2127 else LoadPoints2(a);
2134 void AliTPC::RemapTrackHitIDs(Int_t *map)
2136 if (!fTrackHits) return;
2137 AliObjectArray * arr = fTrackHits->fTrackHitsInfo;
2138 for (UInt_t i=0;i<arr->GetSize();i++){
2139 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2140 info->fTrackID = map[info->fTrackID];
2146 //_____________________________________________________________________________
2147 void AliTPC::LoadPoints2(Int_t)
2150 // Store x, y, z of all hits in memory
2152 if (fTrackHits == 0) return;
2154 Int_t nhits = fTrackHits->GetEntriesFast();
2155 if (nhits == 0) return;
2156 Int_t tracks = gAlice->GetNtrack();
2157 if (fPoints == 0) fPoints = new TObjArray(tracks);
2160 Int_t *ntrk=new Int_t[tracks];
2161 Int_t *limi=new Int_t[tracks];
2162 Float_t **coor=new Float_t*[tracks];
2163 for(Int_t i=0;i<tracks;i++) {
2169 AliPoints *points = 0;
2172 Int_t chunk=nhits/4+1;
2174 // Loop over all the hits and store their position
2176 ahit = FirstHit2(-1);
2177 //for (Int_t hit=0;hit<nhits;hit++) {
2179 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2180 trk=ahit->GetTrack();
2181 if(ntrk[trk]==limi[trk]) {
2183 // Initialise a new track
2184 fp=new Float_t[3*(limi[trk]+chunk)];
2186 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2187 delete [] coor[trk];
2194 fp[3*ntrk[trk] ] = ahit->X();
2195 fp[3*ntrk[trk]+1] = ahit->Y();
2196 fp[3*ntrk[trk]+2] = ahit->Z();
2201 for(trk=0; trk<tracks; ++trk) {
2203 points = new AliPoints();
2204 points->SetMarkerColor(GetMarkerColor());
2205 points->SetMarkerSize(GetMarkerSize());
2206 points->SetDetector(this);
2207 points->SetParticle(trk);
2208 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2209 fPoints->AddAt(points,trk);
2210 delete [] coor[trk];
2220 //_____________________________________________________________________________
2221 void AliTPC::LoadPoints3(Int_t)
2224 // Store x, y, z of all hits in memory
2225 // - only intersection point with pad row
2226 if (fTrackHits == 0) return;
2228 Int_t nhits = fTrackHits->GetEntriesFast();
2229 if (nhits == 0) return;
2230 Int_t tracks = gAlice->GetNtrack();
2231 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2232 fPoints->Expand(2*tracks);
2235 Int_t *ntrk=new Int_t[tracks];
2236 Int_t *limi=new Int_t[tracks];
2237 Float_t **coor=new Float_t*[tracks];
2238 for(Int_t i=0;i<tracks;i++) {
2244 AliPoints *points = 0;
2247 Int_t chunk=nhits/4+1;
2249 // Loop over all the hits and store their position
2251 ahit = FirstHit2(-1);
2252 //for (Int_t hit=0;hit<nhits;hit++) {
2256 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2257 trk=ahit->GetTrack();
2258 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2259 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2260 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2261 if (currentrow!=lastrow){
2262 lastrow = currentrow;
2263 //later calculate intersection point
2264 if(ntrk[trk]==limi[trk]) {
2266 // Initialise a new track
2267 fp=new Float_t[3*(limi[trk]+chunk)];
2269 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2270 delete [] coor[trk];
2277 fp[3*ntrk[trk] ] = ahit->X();
2278 fp[3*ntrk[trk]+1] = ahit->Y();
2279 fp[3*ntrk[trk]+2] = ahit->Z();
2286 for(trk=0; trk<tracks; ++trk) {
2288 points = new AliPoints();
2289 points->SetMarkerColor(GetMarkerColor()+1);
2290 points->SetMarkerStyle(5);
2291 points->SetMarkerSize(0.2);
2292 points->SetDetector(this);
2293 points->SetParticle(trk);
2294 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2295 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2296 fPoints->AddAt(points,tracks+trk);
2297 delete [] coor[trk];
2308 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2312 //fill clones array with intersection of current point with the
2317 const Int_t kcmaxhits=30000;
2318 TVector * xxxx = new TVector(kcmaxhits*4);
2319 TVector & xxx = *xxxx;
2320 Int_t maxhits = kcmaxhits;
2323 AliTPChit * tpcHit=0;
2324 tpcHit = (AliTPChit*)FirstHit2(-1);
2325 Int_t currentIndex=0;
2326 Int_t lastrow=-1; //last writen row
2329 if (tpcHit==0) continue;
2330 sector=tpcHit->fSector; // sector number
2331 ipart=tpcHit->Track();
2335 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2336 Int_t index[3]={1,sector,0};
2337 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2338 if (currentrow<0) continue;
2339 if (lastrow<0) lastrow=currentrow;
2340 if (currentrow==lastrow){
2341 if ( currentIndex>=maxhits){
2343 xxx.ResizeTo(4*maxhits);
2345 xxx(currentIndex*4)=x[0];
2346 xxx(currentIndex*4+1)=x[1];
2347 xxx(currentIndex*4+2)=x[2];
2348 xxx(currentIndex*4+3)=tpcHit->fQ;
2352 if (currentIndex>2){
2364 for (Int_t index=0;index<currentIndex;index++){
2366 x=x2=x3=x4=xxx(index*4);
2374 sumy+=xxx(index*4+1);
2375 sumxy+=xxx(index*4+1)*x;
2376 sumx2y+=xxx(index*4+1)*x2;
2377 sumz+=xxx(index*4+2);
2378 sumxz+=xxx(index*4+2)*x;
2379 sumx2z+=xxx(index*4+2)*x2;
2380 sumq+=xxx(index*4+3);
2382 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2383 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2384 sumx2*(sumx*sumx3-sumx2*sumx2);
2386 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2387 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2388 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2389 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2391 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2392 sumx2*(sumx*sumx2y-sumx2*sumxy);
2393 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2394 sumx2*(sumx*sumx2z-sumx2*sumxz);
2396 Float_t y=detay/det+centralPad;
2397 Float_t z=detaz/det;
2398 Float_t by=detby/det; //y angle
2399 Float_t bz=detbz/det; //z angle
2400 sumy/=Float_t(currentIndex);
2401 sumz/=Float_t(currentIndex);
2403 AliComplexCluster cl;
2409 cl.fTracks[0]=ipart;
2411 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2412 if (row!=0) row->InsertCluster(&cl);
2415 } //end of calculating cluster for given row
2417 } // end of loop over hits