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.37 2001/06/12 07:17:18 kowal2
19 Hits2SDigits method implemented (summable digits)
21 Revision 1.36 2001/05/16 14:57:25 alibrary
22 New files for folders and Stack
24 Revision 1.35 2001/05/08 16:02:22 kowal2
25 Updated material specifications
27 Revision 1.34 2001/05/08 15:00:15 hristov
28 Corrections for tracking in arbitrary magnenetic field. Changes towards a concept of global Alice track. Back propagation of reconstructed tracks (Yu.Belikov)
30 Revision 1.33 2001/04/03 12:40:43 kowal2
33 Revision 1.32 2001/03/12 17:47:36 hristov
34 Changes needed on Sun with CC 5.0
36 Revision 1.31 2001/03/12 08:21:50 kowal2
37 Corrected C++ bug in the material definitions
39 Revision 1.30 2001/03/01 17:34:47 kowal2
40 Correction due to the accuracy problem
42 Revision 1.29 2001/02/28 16:34:40 kowal2
43 Protection against nonphysical values of the avalanche size,
46 Revision 1.28 2001/01/26 19:57:19 hristov
47 Major upgrade of AliRoot code
49 Revision 1.27 2001/01/13 17:29:33 kowal2
50 Sun compiler correction
52 Revision 1.26 2001/01/10 07:59:43 kowal2
53 Corrections to load points from the noncompressed hits.
55 Revision 1.25 2000/11/02 07:25:31 kowal2
56 Changes due to the new hit structure.
59 Revision 1.24 2000/10/05 16:06:09 kowal2
60 Forward declarations. Changes due to a new class AliComplexCluster.
62 Revision 1.23 2000/10/02 21:28:18 fca
63 Removal of useless dependecies via forward declarations
65 Revision 1.22 2000/07/10 20:57:39 hristov
66 Update of TPC code and macros by M.Kowalski
68 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
69 Changes to obey the coding rules
71 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
72 Splitted from AliTPCtracking
74 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
75 Changed parameter settings
77 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
79 Defaults loaded automatically (hard-wired)
80 Optional parameters can be set via macro called in the constructor
82 Revision 1.19 2000/04/18 19:00:59 fca
83 Small bug fixes to TPC files
85 Revision 1.18 2000/04/17 09:37:33 kowal2
86 removed obsolete AliTPCDigitsDisplay.C
88 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
90 New, experimental data structure from M. Ivanov
91 New tracking algorithm
92 Different pad geometry for different sectors
93 Digitization rewritten
95 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
96 Not used anymore - removed
98 Revision 1.17 2000/01/19 17:17:30 fca
99 Introducing a list of lists of hits -- more hits allowed for detector now
101 Revision 1.16 1999/11/05 09:29:23 fca
102 Accept only signals > 0
104 Revision 1.15 1999/10/08 06:26:53 fca
105 Removed ClustersIndex - not used anymore
107 Revision 1.14 1999/09/29 09:24:33 fca
108 Introduction of the Copyright and cvs Log
112 ///////////////////////////////////////////////////////////////////////////////
114 // Time Projection Chamber //
115 // This class contains the basic functions for the Time Projection Chamber //
116 // detector. Functions specific to one particular geometry are //
117 // contained in the derived classes //
121 <img src="picts/AliTPCClass.gif">
126 ///////////////////////////////////////////////////////////////////////////////
134 #include <TGeometry.h>
137 #include <TObjectTable.h>
138 #include "TParticle.h"
144 #include <iostream.h>
151 #include "AliTPCParamSR.h"
152 #include "AliTPCPRF2D.h"
153 #include "AliTPCRF1D.h"
154 #include "AliDigits.h"
155 #include "AliSimDigits.h"
156 #include "AliTPCTrackHits.h"
157 #include "AliPoints.h"
158 #include "AliArrayBranch.h"
161 #include "AliTPCDigitsArray.h"
162 #include "AliComplexCluster.h"
163 #include "AliClusters.h"
164 #include "AliTPCClustersRow.h"
165 #include "AliTPCClustersArray.h"
167 #include "AliTPCcluster.h"
168 #include "AliTPCclusterer.h"
169 #include "AliTPCtracker.h"
171 #include <TInterpreter.h>
178 //_____________________________________________________________________________
182 // Default constructor
197 //_____________________________________________________________________________
198 AliTPC::AliTPC(const char *name, const char *title)
199 : AliDetector(name,title)
202 // Standard constructor
206 // Initialise arrays of hits and digits
207 fHits = new TClonesArray("AliTPChit", 176);
208 gAlice->AddHitList(fHits);
214 fTrackHits = new AliTPCTrackHits; //MI - 13.09.2000
215 fTrackHits->SetHitPrecision(0.002);
216 fTrackHits->SetStepPrecision(0.003);
217 fTrackHits->SetMaxDistance(100);
221 // Initialise counters
227 // Initialise color attributes
228 SetMarkerColor(kYellow);
231 // Set TPC parameters
235 if (!strcmp(title,"Default")) {
236 fTPCParam = new AliTPCParamSR;
238 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
244 //_____________________________________________________________________________
255 delete fTrackHits; //MI 15.09.2000
258 //_____________________________________________________________________________
259 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
262 // Add a hit to the list
264 // TClonesArray &lhits = *fHits;
265 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
267 TClonesArray &lhits = *fHits;
268 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
271 AddHit2(track,vol,hits);
274 //_____________________________________________________________________________
275 void AliTPC::BuildGeometry()
279 // Build TPC ROOT TNode geometry for the event display
284 const int kColorTPC=19;
285 char name[5], title[25];
286 const Double_t kDegrad=TMath::Pi()/180;
287 const Double_t kRaddeg=180./TMath::Pi();
290 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
291 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
293 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
294 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
296 Int_t nLo = fTPCParam->GetNInnerSector()/2;
297 Int_t nHi = fTPCParam->GetNOuterSector()/2;
299 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
300 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
301 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
302 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
305 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
306 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
312 // Get ALICE top node
315 nTop=gAlice->GetGeometry()->GetNode("alice");
319 rl = fTPCParam->GetInnerRadiusLow();
320 ru = fTPCParam->GetInnerRadiusUp();
324 sprintf(name,"LS%2.2d",i);
326 sprintf(title,"TPC low sector %3d",i);
329 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
330 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
331 tubs->SetNumberOfDivisions(1);
333 nNode = new TNode(name,title,name,0,0,0,"");
334 nNode->SetLineColor(kColorTPC);
340 rl = fTPCParam->GetOuterRadiusLow();
341 ru = fTPCParam->GetOuterRadiusUp();
344 sprintf(name,"US%2.2d",i);
346 sprintf(title,"TPC upper sector %d",i);
348 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
349 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
350 tubs->SetNumberOfDivisions(1);
352 nNode = new TNode(name,title,name,0,0,0,"");
353 nNode->SetLineColor(kColorTPC);
359 //_____________________________________________________________________________
360 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
363 // Calculate distance from TPC to mouse on the display
369 void AliTPC::Clusters2Tracks(TFile *of) {
370 //-----------------------------------------------------------------
371 // This is a track finder.
372 //-----------------------------------------------------------------
373 AliTPCtracker tracker(fTPCParam);
374 tracker.Clusters2Tracks(gFile,of);
377 //_____________________________________________________________________________
378 void AliTPC::CreateMaterials()
380 //-----------------------------------------------
381 // Create Materials for for TPC simulations
382 //-----------------------------------------------
384 //-----------------------------------------------------------------
385 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
386 //-----------------------------------------------------------------
388 Int_t iSXFLD=gAlice->Field()->Integ();
389 Float_t sXMGMX=gAlice->Field()->Max();
391 Float_t amat[5]; // atomic numbers
392 Float_t zmat[5]; // z
393 Float_t wmat[5]; // proportions
399 //***************** Gases *************************
401 //-------------------------------------------------
403 //-------------------------------------------------
414 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
424 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
427 //--------------------------------------------------------------
429 //--------------------------------------------------------------
446 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
448 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
463 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
465 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
481 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
483 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
485 //----------------------------------------------------------------
486 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
487 //----------------------------------------------------------------
493 Float_t rho,absl,X0,buf[1];
497 for(nc = 0;nc<fNoComp;nc++)
500 // retrive material constants
502 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
507 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
509 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
510 density += fMixtProp[nc]*rho; // density of the mixture
514 // mixture proportions by weight!
516 for(nc = 0;nc<fNoComp;nc++)
519 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
521 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
522 apure[nnc] : amol[nnc])/am;
526 // Drift gases 1 - nonsensitive, 2 - sensitive
528 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
529 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
538 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
541 //----------------------------------------------------------------------
543 //----------------------------------------------------------------------
565 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
587 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
605 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
623 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
625 // SiO2 - used later for the glass fiber
637 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
646 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
655 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
664 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
682 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
701 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
703 // Epoxy - C14 H20 O3
720 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
728 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
732 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,X0,absl,buf,nbuf);
736 wmat[0]=0.644; // by weight!
739 density=0.5*(1.25+2.265);
741 AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
745 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
747 wmat[0]=0.725; // by weight!
752 AliMixture(39,"G10",amat,zmat,density,2,wmat);
757 //----------------------------------------------------------
758 // tracking media for gases
759 //----------------------------------------------------------
761 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
762 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
763 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
764 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
766 //-----------------------------------------------------------
767 // tracking media for solids
768 //-----------------------------------------------------------
770 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
771 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
772 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
773 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
774 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
775 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
776 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
777 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
778 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
779 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
780 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
781 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
786 void AliTPC::Digits2Clusters(TFile *of, Int_t eventnumber)
788 //-----------------------------------------------------------------
789 // This is a simple cluster finder.
790 //-----------------------------------------------------------------
791 AliTPCclusterer::Digits2Clusters(fTPCParam,of,eventnumber);
794 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
795 extern Double_t SigmaZ2(Double_t, Double_t);
796 //_____________________________________________________________________________
797 void AliTPC::Hits2Clusters(TFile *of, Int_t eventn)
799 //--------------------------------------------------------
800 // TPC simple cluster generator from hits
801 // obtained from the TPC Fast Simulator
802 // The point errors are taken from the parametrization
803 //--------------------------------------------------------
805 //-----------------------------------------------------------------
806 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
807 //-----------------------------------------------------------------
808 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
809 // Jouri.Belikov@cern.ch
810 //----------------------------------------------------------------
812 /////////////////////////////////////////////////////////////////////////////
814 //---------------------------------------------------------------------
815 // ALICE TPC Cluster Parameters
816 //--------------------------------------------------------------------
820 // Cluster width in rphi
821 const Float_t kACrphi=0.18322;
822 const Float_t kBCrphi=0.59551e-3;
823 const Float_t kCCrphi=0.60952e-1;
824 // Cluster width in z
825 const Float_t kACz=0.19081;
826 const Float_t kBCz=0.55938e-3;
827 const Float_t kCCz=0.30428;
829 TDirectory *savedir=gDirectory;
832 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
836 //if(fDefaults == 0) SetDefaults();
838 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
840 TParticle *particle; // pointer to a given particle
841 AliTPChit *tpcHit; // pointer to a sigle TPC hit
845 Float_t pl,pt,tanth,rpad,ratio;
848 //---------------------------------------------------------------
849 // Get the access to the tracks
850 //---------------------------------------------------------------
852 TTree *tH = gAlice->TreeH();
853 Stat_t ntracks = tH->GetEntries();
855 //Switch to the output file
860 sprintf(cname,"TreeC_TPC_%d",eventn);
862 fTPCParam->Write(fTPCParam->GetTitle());
863 AliTPCClustersArray carray;
864 carray.Setup(fTPCParam);
865 carray.SetClusterType("AliTPCcluster");
868 Int_t nclusters=0; //cluster counter
870 //------------------------------------------------------------
871 // Loop over all sectors (72 sectors for 20 deg
872 // segmentation for both lower and upper sectors)
873 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
874 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
876 // First cluster for sector 0 starts at "0"
877 //------------------------------------------------------------
879 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
881 fTPCParam->AdjustCosSin(isec,cph,sph);
883 //------------------------------------------------------------
885 //------------------------------------------------------------
887 for(Int_t track=0;track<ntracks;track++){
891 // Get number of the TPC hits
893 // nhits=fHits->GetEntriesFast();
896 tpcHit = (AliTPChit*)FirstHit(-1);
900 // for(Int_t hit=0;hit<nhits;hit++){
901 //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
905 if (tpcHit->fQ == 0.) {
906 tpcHit = (AliTPChit*) NextHit();
907 continue; //information about track (I.Belikov)
909 sector=tpcHit->fSector; // sector number
912 // if(sector != isec) continue; //terminate iteration
915 tpcHit = (AliTPChit*) NextHit();
918 ipart=tpcHit->Track();
919 particle=gAlice->Particle(ipart);
922 if(pt < 1.e-9) pt=1.e-9;
924 tanth = TMath::Abs(tanth);
925 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
926 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
928 // space-point resolutions
930 sigmaRphi=SigmaY2(rpad,tanth,pt);
931 sigmaZ =SigmaZ2(rpad,tanth );
935 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
936 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
938 // temporary protection
940 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
941 if(sigmaZ < 0.) sigmaZ=0.4e-3;
942 if(clRphi < 0.) clRphi=2.5e-3;
943 if(clZ < 0.) clZ=2.5e-5;
948 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
949 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
951 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
952 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
953 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
954 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
955 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
956 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
957 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
958 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
959 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
960 xyz[2]=sigmaRphi; // fSigmaY2
961 xyz[3]=sigmaZ; // fSigmaZ2
962 xyz[4]=tpcHit->fQ; // q
964 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
965 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
967 Int_t tracks[3]={tpcHit->Track(), -1, -1};
968 AliTPCcluster cluster(tracks,xyz);
970 clrow->InsertCluster(&cluster); nclusters++;
972 tpcHit = (AliTPChit*)NextHit();
975 } // end of loop over hits
977 } // end of loop over tracks
979 Int_t nrows=fTPCParam->GetNRow(isec);
980 for (Int_t irow=0; irow<nrows; irow++) {
981 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
982 if (!clrow) continue;
983 carray.StoreRow(isec,irow);
984 carray.ClearRow(isec,irow);
987 } // end of loop over sectors
989 cerr<<"Number of made clusters : "<<nclusters<<" \n";
991 carray.GetTree()->Write(cname);
993 savedir->cd(); //switch back to the input file
997 //_________________________________________________________________
998 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1000 //--------------------------------------------------------
1001 //calculate exact cross point of track and given pad row
1002 //resulting values are expressed in "digit" coordinata
1003 //--------------------------------------------------------
1005 //-----------------------------------------------------------------
1006 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1007 //-----------------------------------------------------------------
1009 if (fClustersArray==0){
1013 TParticle *particle; // pointer to a given particle
1014 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1017 const Int_t kcmaxhits=30000;
1018 TVector * xxxx = new TVector(kcmaxhits*4);
1019 TVector & xxx = *xxxx;
1020 Int_t maxhits = kcmaxhits;
1021 //construct array for each padrow
1022 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
1023 fClustersArray->CreateRow(isec,i);
1025 //---------------------------------------------------------------
1026 // Get the access to the tracks
1027 //---------------------------------------------------------------
1029 TTree *tH = gAlice->TreeH();
1030 Stat_t ntracks = tH->GetEntries();
1031 Int_t npart = gAlice->GetNtrack();
1033 //------------------------------------------------------------
1035 //------------------------------------------------------------
1037 for(Int_t track=0;track<ntracks;track++){
1039 tH->GetEvent(track);
1041 // Get number of the TPC hits and a pointer
1044 nhits=fHits->GetEntriesFast();
1048 Int_t currentIndex=0;
1049 Int_t lastrow=-1; //last writen row
1050 for(Int_t hit=0;hit<nhits;hit++){
1051 tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
1052 if (tpcHit==0) continue;
1053 sector=tpcHit->fSector; // sector number
1054 if(sector != isec) continue;
1055 ipart=tpcHit->Track();
1056 if (ipart<npart) particle=gAlice->Particle(ipart);
1060 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1061 Int_t index[3]={1,isec,0};
1062 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1063 if (currentrow<0) continue;
1064 if (lastrow<0) lastrow=currentrow;
1065 if (currentrow==lastrow){
1066 if ( currentIndex>=maxhits){
1068 xxx.ResizeTo(4*maxhits);
1070 xxx(currentIndex*4)=x[0];
1071 xxx(currentIndex*4+1)=x[1];
1072 xxx(currentIndex*4+2)=x[2];
1073 xxx(currentIndex*4+3)=tpcHit->fQ;
1077 if (currentIndex>2){
1089 for (Int_t index=0;index<currentIndex;index++){
1091 x=x2=x3=x4=xxx(index*4);
1099 sumy+=xxx(index*4+1);
1100 sumxy+=xxx(index*4+1)*x;
1101 sumx2y+=xxx(index*4+1)*x2;
1102 sumz+=xxx(index*4+2);
1103 sumxz+=xxx(index*4+2)*x;
1104 sumx2z+=xxx(index*4+2)*x2;
1105 sumq+=xxx(index*4+3);
1107 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1108 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1109 sumx2*(sumx*sumx3-sumx2*sumx2);
1111 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1112 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1113 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1114 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1116 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1117 sumx2*(sumx*sumx2y-sumx2*sumxy);
1118 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1119 sumx2*(sumx*sumx2z-sumx2*sumxz);
1121 Float_t y=detay/det+centralPad;
1122 Float_t z=detaz/det;
1123 Float_t by=detby/det; //y angle
1124 Float_t bz=detbz/det; //z angle
1125 sumy/=Float_t(currentIndex);
1126 sumz/=Float_t(currentIndex);
1127 AliComplexCluster cl;
1133 cl.fTracks[0]=ipart;
1135 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1136 if (row!=0) row->InsertCluster(&cl);
1139 } //end of calculating cluster for given row
1143 } // end of loop over hits
1144 } // end of loop over tracks
1145 //write padrows to tree
1146 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1147 fClustersArray->StoreRow(isec,ii);
1148 fClustersArray->ClearRow(isec,ii);
1153 //___________________________________________
1154 void AliTPC::SDigits2Digits(Int_t eventnumber)
1158 cerr<<"Digitizing TPC...\n";
1160 Hits2Digits(eventnumber);
1165 // char treeName[100];
1167 // sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1169 // GetDigitsArray()->GetTree()->Write(treeName);
1171 //__________________________________________________________________
1172 void AliTPC::SetDefaults(){
1175 cerr<<"Setting default parameters...\n";
1177 // Set response functions
1179 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1180 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1181 AliTPCPRF2D * prfouter = new AliTPCPRF2D;
1182 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1183 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1184 rf->SetOffset(3*param->GetZSigma());
1187 TDirectory *savedir=gDirectory;
1188 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1190 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1193 prfinner->Read("prf_07504_Gati_056068_d02");
1194 prfouter->Read("prf_10006_Gati_047051_d03");
1198 param->SetInnerPRF(prfinner);
1199 param->SetOuterPRF(prfouter);
1200 param->SetTimeRF(rf);
1210 //__________________________________________________________________
1211 void AliTPC::Hits2Digits(Int_t eventnumber)
1213 //----------------------------------------------------
1214 // Loop over all sectors for a single event
1215 //----------------------------------------------------
1218 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1220 //setup TPCDigitsArray
1222 if(GetDigitsArray()) delete GetDigitsArray();
1224 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1225 arr->SetClass("AliSimDigits");
1226 arr->Setup(fTPCParam);
1227 arr->MakeTree(fDigitsFile);
1228 SetDigitsArray(arr);
1230 fDigitsSwitch=0; // standard digits
1232 cerr<<"Digitizing TPC...\n";
1234 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1240 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1242 GetDigitsArray()->GetTree()->Write(treeName);
1249 //__________________________________________________________________
1250 void AliTPC::Hits2SDigits(Int_t eventnumber)
1253 //-----------------------------------------------------------
1254 // summable digits - 16 bit "ADC", no noise, no saturation
1255 //-----------------------------------------------------------
1257 //----------------------------------------------------
1258 // Loop over all sectors for a single event
1259 //----------------------------------------------------
1262 if(fDefaults == 0) SetDefaults();
1264 //setup TPCDigitsArray
1266 if(GetDigitsArray()) delete GetDigitsArray();
1268 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1269 arr->SetClass("AliSimDigits");
1270 arr->Setup(fTPCParam);
1271 arr->MakeTree(fDigitsFile);
1272 SetDigitsArray(arr);
1274 cerr<<"Digitizing TPC...\n";
1276 fDigitsSwitch=1; // summable digits
1278 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1285 sprintf(treeName,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1287 GetDigitsArray()->GetTree()->Write(treeName);
1293 //_____________________________________________________________________________
1294 void AliTPC::Hits2DigitsSector(Int_t isec)
1296 //-------------------------------------------------------------------
1297 // TPC conversion from hits to digits.
1298 //-------------------------------------------------------------------
1300 //-----------------------------------------------------------------
1301 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1302 //-----------------------------------------------------------------
1304 //-------------------------------------------------------
1305 // Get the access to the track hits
1306 //-------------------------------------------------------
1308 // check if the parameters are set - important if one calls this method
1309 // directly, not from the Hits2Digits
1311 if(fDefaults == 0) SetDefaults();
1313 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1314 Stat_t ntracks = tH->GetEntries();
1318 //-------------------------------------------
1319 // Only if there are any tracks...
1320 //-------------------------------------------
1324 printf("*** Processing sector number %d ***\n",isec);
1326 Int_t nrows =fTPCParam->GetNRow(isec);
1328 row= new TObjArray* [nrows];
1330 MakeSector(isec,nrows,tH,ntracks,row);
1332 //--------------------------------------------------------
1333 // Digitize this sector, row by row
1334 // row[i] is the pointer to the TObjArray of TVectors,
1335 // each one containing electrons accepted on this
1336 // row, assigned into tracks
1337 //--------------------------------------------------------
1341 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1343 for (i=0;i<nrows;i++){
1345 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1347 DigitizeRow(i,isec,row);
1349 fDigitsArray->StoreRow(isec,i);
1351 Int_t ndig = dig->GetDigitSize();
1354 printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1356 fDigitsArray->ClearRow(isec,i);
1359 } // end of the sector digitization
1361 for(i=0;i<nrows;i++){
1366 delete [] row; // delete the array of pointers to TObjArray-s
1370 } // end of Hits2DigitsSector
1373 //_____________________________________________________________________________
1374 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1376 //-----------------------------------------------------------
1377 // Single row digitization, coupling from the neighbouring
1378 // rows taken into account
1379 //-----------------------------------------------------------
1381 //-----------------------------------------------------------------
1382 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1383 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1384 //-----------------------------------------------------------------
1387 Float_t zerosup = fTPCParam->GetZeroSup();
1388 Int_t nrows =fTPCParam->GetNRow(isec);
1389 fCurrentIndex[1]= isec;
1392 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1393 Int_t nofTbins = fTPCParam->GetMaxTBin();
1394 Int_t indexRange[4];
1396 // Integrated signal for this row
1397 // and a single track signal
1399 TMatrix *m1 = new TMatrix(0,nofPads,0,nofTbins); // integrated
1400 TMatrix *m2 = new TMatrix(0,nofPads,0,nofTbins); // single
1402 TMatrix &total = *m1;
1404 // Array of pointers to the label-signal list
1406 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1407 Float_t **pList = new Float_t* [nofDigits];
1411 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1415 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1416 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1417 for (Int_t row= row1;row<=row2;row++){
1418 Int_t nTracks= rows[row]->GetEntries();
1419 for (i1=0;i1<nTracks;i1++){
1420 fCurrentIndex[2]= row;
1421 fCurrentIndex[3]=irow;
1423 m2->Zero(); // clear single track signal matrix
1424 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1425 GetList(trackLabel,nofPads,m2,indexRange,pList);
1427 else GetSignal(rows[row],i1,0,m1,indexRange);
1433 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1434 for(Int_t ip=0;ip<nofPads;ip++){
1435 for(Int_t it=0;it<nofTbins;it++){
1437 Float_t q = total(ip,it);
1439 Int_t gi =it*nofPads+ip; // global index
1441 if(fDigitsSwitch == 0){
1443 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
1447 if(q <=zerosup) continue; // do not fill zeros
1448 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1458 // "real" signal or electronic noise (list = -1)?
1461 for(Int_t j1=0;j1<3;j1++){
1462 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1467 <A NAME="AliDigits"></A>
1468 using of AliDigits object
1471 dig->SetDigitFast((Short_t)q,it,ip);
1472 if (fDigitsArray->IsSimulated())
1474 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1475 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1476 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1480 } // end of loop over time buckets
1481 } // end of lop over pads
1484 // This row has been digitized, delete nonused stuff
1487 for(lp=0;lp<nofDigits;lp++){
1488 if(pList[lp]) delete [] pList[lp];
1497 } // end of DigitizeRow
1499 //_____________________________________________________________________________
1501 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, TMatrix *m1, TMatrix *m2,
1505 //---------------------------------------------------------------
1506 // Calculates 2-D signal (pad,time) for a single track,
1507 // returns a pointer to the signal matrix and the track label
1508 // No digitization is performed at this level!!!
1509 //---------------------------------------------------------------
1511 //-----------------------------------------------------------------
1512 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1513 // Modified: Marian Ivanov
1514 //-----------------------------------------------------------------
1518 tv = (TVector*)p1->At(ntr); // pointer to a track
1521 Float_t label = v(0);
1522 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1524 Int_t nElectrons = (tv->GetNrows()-1)/4;
1525 indexRange[0]=9999; // min pad
1526 indexRange[1]=-1; // max pad
1527 indexRange[2]=9999; //min time
1528 indexRange[3]=-1; // max time
1530 // Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1532 TMatrix &signal = *m1;
1533 TMatrix &total = *m2;
1535 // Loop over all electrons
1537 for(Int_t nel=0; nel<nElectrons; nel++){
1539 Float_t aval = v(idx+4);
1540 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1541 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1542 Int_t n = fTPCParam->CalcResponse(xyz,fCurrentIndex,fCurrentIndex[3]);
1544 if (n>0) for (Int_t i =0; i<n; i++){
1545 Int_t *index = fTPCParam->GetResBin(i);
1546 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1547 if ( ( pad<(fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]))) && (pad>0)) {
1548 Int_t time=index[2];
1549 Float_t weight = fTPCParam->GetResWeight(i); //we normalise response to ADC channel
1550 weight *= eltoadcfac;
1552 if (m1!=0) signal(pad,time)+=weight;
1553 total(pad,time)+=weight;
1554 indexRange[0]=TMath::Min(indexRange[0],pad);
1555 indexRange[1]=TMath::Max(indexRange[1],pad);
1556 indexRange[2]=TMath::Min(indexRange[2],time);
1557 indexRange[3]=TMath::Max(indexRange[3],time);
1560 } // end of loop over electrons
1562 return label; // returns track label when finished
1565 //_____________________________________________________________________________
1566 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *indexRange,
1569 //----------------------------------------------------------------------
1570 // Updates the list of tracks contributing to digits for a given row
1571 //----------------------------------------------------------------------
1573 //-----------------------------------------------------------------
1574 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1575 //-----------------------------------------------------------------
1577 TMatrix &signal = *m;
1579 // lop over nonzero digits
1581 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1582 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1585 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1587 if(signal(ip,it)<0.5) continue;
1590 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1592 if(!pList[globalIndex]){
1595 // Create new list (6 elements - 3 signals and 3 labels),
1598 pList[globalIndex] = new Float_t [6];
1602 *pList[globalIndex] = -1.;
1603 *(pList[globalIndex]+1) = -1.;
1604 *(pList[globalIndex]+2) = -1.;
1605 *(pList[globalIndex]+3) = -1.;
1606 *(pList[globalIndex]+4) = -1.;
1607 *(pList[globalIndex]+5) = -1.;
1610 *pList[globalIndex] = label;
1611 *(pList[globalIndex]+3) = signal(ip,it);
1615 // check the signal magnitude
1617 Float_t highest = *(pList[globalIndex]+3);
1618 Float_t middle = *(pList[globalIndex]+4);
1619 Float_t lowest = *(pList[globalIndex]+5);
1622 // compare the new signal with already existing list
1625 if(signal(ip,it)<lowest) continue; // neglect this track
1629 if (signal(ip,it)>highest){
1630 *(pList[globalIndex]+5) = middle;
1631 *(pList[globalIndex]+4) = highest;
1632 *(pList[globalIndex]+3) = signal(ip,it);
1634 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1635 *(pList[globalIndex]+1) = *pList[globalIndex];
1636 *pList[globalIndex] = label;
1638 else if (signal(ip,it)>middle){
1639 *(pList[globalIndex]+5) = middle;
1640 *(pList[globalIndex]+4) = signal(ip,it);
1642 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1643 *(pList[globalIndex]+1) = label;
1646 *(pList[globalIndex]+5) = signal(ip,it);
1647 *(pList[globalIndex]+2) = label;
1651 } // end of loop over pads
1652 } // end of loop over time bins
1657 //___________________________________________________________________
1658 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1659 Stat_t ntracks,TObjArray **row)
1662 //-----------------------------------------------------------------
1663 // Prepares the sector digitization, creates the vectors of
1664 // tracks for each row of this sector. The track vector
1665 // contains the track label and the position of electrons.
1666 //-----------------------------------------------------------------
1668 //-----------------------------------------------------------------
1669 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1670 //-----------------------------------------------------------------
1672 Float_t gasgain = fTPCParam->GetGasGain();
1676 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1679 if (fHitType&2) branch = TH->GetBranch("TPC2");
1680 else branch = TH->GetBranch("TPC");
1683 //----------------------------------------------
1684 // Create TObjArray-s, one for each row,
1685 // each TObjArray will store the TVectors
1686 // of electrons, one TVector per each track.
1687 //----------------------------------------------
1689 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
1690 TVector **tracks = new TVector* [nrows]; //pointers to the track vectors
1691 for(i=0; i<nrows; i++){
1692 row[i] = new TObjArray;
1699 //--------------------------------------------------------------------
1700 // Loop over tracks, the "track" contains the full history
1701 //--------------------------------------------------------------------
1703 Int_t previousTrack,currentTrack;
1704 previousTrack = -1; // nothing to store so far!
1706 for(Int_t track=0;track<ntracks;track++){
1707 Bool_t isInSector=kTRUE;
1712 TBranch * br = TH->GetBranch("fTrackHitsInfo");
1713 br->GetEvent(track);
1714 AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
1715 for (UInt_t j=0;j<ar->GetSize();j++){
1716 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
1719 if (!isInSector) continue;
1721 branch->GetEntry(track); // get next track
1725 tpcHit = (AliTPChit*)FirstHit(-1);
1727 //--------------------------------------------------------------
1729 //--------------------------------------------------------------
1734 Int_t sector=tpcHit->fSector; // sector number
1735 // if(sector != isec) continue;
1737 tpcHit = (AliTPChit*) NextHit();
1741 currentTrack = tpcHit->Track(); // track number
1744 if(currentTrack != previousTrack){
1746 // store already filled fTrack
1748 for(i=0;i<nrows;i++){
1749 if(previousTrack != -1){
1750 if(nofElectrons[i]>0){
1751 TVector &v = *tracks[i];
1752 v(0) = previousTrack;
1753 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1754 row[i]->Add(tracks[i]);
1757 delete tracks[i]; // delete empty TVector
1763 tracks[i] = new TVector(481); // TVectors for the next fTrack
1765 } // end of loop over rows
1767 previousTrack=currentTrack; // update track label
1770 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1772 //---------------------------------------------------
1773 // Calculate the electron attachment probability
1774 //---------------------------------------------------
1777 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
1778 /fTPCParam->GetDriftV();
1780 Float_t attProb = fTPCParam->GetAttCoef()*
1781 fTPCParam->GetOxyCont()*time; // fraction!
1783 //-----------------------------------------------
1784 // Loop over electrons
1785 //-----------------------------------------------
1788 for(Int_t nel=0;nel<qI;nel++){
1789 // skip if electron lost due to the attachment
1790 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1795 // protection for the nonphysical avalanche size (10**6 maximum)
1797 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
1798 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
1801 TransportElectron(xyz,index); //MI change -august
1803 fTPCParam->GetPadRow(xyz,index); //MI change august
1804 rowNumber = index[2];
1805 //transform position to local digit coordinates
1806 //relative to nearest pad row
1807 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
1808 nofElectrons[rowNumber]++;
1809 //----------------------------------
1810 // Expand vector if necessary
1811 //----------------------------------
1812 if(nofElectrons[rowNumber]>120){
1813 Int_t range = tracks[rowNumber]->GetNrows();
1814 if((nofElectrons[rowNumber])>(range-1)/4){
1816 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
1820 TVector &v = *tracks[rowNumber];
1821 Int_t idx = 4*nofElectrons[rowNumber]-3;
1823 v(idx)= xyz[0]; // X - pad row coordinate
1824 v(idx+1)=xyz[1]; // Y - pad coordinate (along the pad-row)
1825 v(idx+2)=xyz[2]; // Z - time bin coordinate
1826 v(idx+3)=xyz[3]; // avalanche size
1827 } // end of loop over electrons
1829 tpcHit = (AliTPChit*)NextHit();
1831 } // end of loop over hits
1832 } // end of loop over tracks
1835 // store remaining track (the last one) if not empty
1838 for(i=0;i<nrows;i++){
1839 if(nofElectrons[i]>0){
1840 TVector &v = *tracks[i];
1841 v(0) = previousTrack;
1842 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1843 row[i]->Add(tracks[i]);
1852 delete [] nofElectrons;
1855 } // end of MakeSector
1858 //_____________________________________________________________________________
1862 // Initialise TPC detector after definition of geometry
1867 printf("\n%s: ",ClassName());
1868 for(i=0;i<35;i++) printf("*");
1869 printf(" TPC_INIT ");
1870 for(i=0;i<35;i++) printf("*");
1871 printf("\n%s: ",ClassName());
1873 for(i=0;i<80;i++) printf("*");
1878 //_____________________________________________________________________________
1879 void AliTPC::MakeBranch(Option_t* option, const char *file)
1882 // Create Tree branches for the TPC.
1884 Int_t buffersize = 4000;
1885 char branchname[10];
1886 sprintf(branchname,"%s",GetName());
1888 AliDetector::MakeBranch(option,file);
1890 const char *d = strstr(option,"D");
1892 if (fDigits && gAlice->TreeD() && d) {
1893 MakeBranchInTree(gAlice->TreeD(),
1894 branchname, &fDigits, buffersize, file);
1897 if (fHitType&2) MakeBranch2(option,file); // MI change 14.09.2000
1900 //_____________________________________________________________________________
1901 void AliTPC::ResetDigits()
1904 // Reset number of digits and the digits array for this detector
1907 if (fDigits) fDigits->Clear();
1910 //_____________________________________________________________________________
1911 void AliTPC::SetSecAL(Int_t sec)
1913 //---------------------------------------------------
1914 // Activate/deactivate selection for lower sectors
1915 //---------------------------------------------------
1917 //-----------------------------------------------------------------
1918 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1919 //-----------------------------------------------------------------
1924 //_____________________________________________________________________________
1925 void AliTPC::SetSecAU(Int_t sec)
1927 //----------------------------------------------------
1928 // Activate/deactivate selection for upper sectors
1929 //---------------------------------------------------
1931 //-----------------------------------------------------------------
1932 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1933 //-----------------------------------------------------------------
1938 //_____________________________________________________________________________
1939 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
1941 //----------------------------------------
1942 // Select active lower sectors
1943 //----------------------------------------
1945 //-----------------------------------------------------------------
1946 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1947 //-----------------------------------------------------------------
1957 //_____________________________________________________________________________
1958 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
1959 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
1960 Int_t s11 , Int_t s12)
1962 //--------------------------------
1963 // Select active upper sectors
1964 //--------------------------------
1966 //-----------------------------------------------------------------
1967 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1968 //-----------------------------------------------------------------
1984 //_____________________________________________________________________________
1985 void AliTPC::SetSens(Int_t sens)
1988 //-------------------------------------------------------------
1989 // Activates/deactivates the sensitive strips at the center of
1990 // the pad row -- this is for the space-point resolution calculations
1991 //-------------------------------------------------------------
1993 //-----------------------------------------------------------------
1994 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1995 //-----------------------------------------------------------------
2001 void AliTPC::SetSide(Float_t side=0.)
2003 // choice of the TPC side
2008 //____________________________________________________________________________
2009 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2010 Float_t p2,Float_t p3)
2013 // gax mixture definition
2027 //_____________________________________________________________________________
2029 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2032 // electron transport taking into account:
2034 // 2.ExB at the wires
2035 // 3. nonisochronity
2037 // xyz and index must be already transformed to system 1
2040 fTPCParam->Transform1to2(xyz,index);
2043 Float_t driftl=xyz[2];
2044 if(driftl<0.01) driftl=0.01;
2045 driftl=TMath::Sqrt(driftl);
2046 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2047 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2048 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2049 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2050 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2054 if (fTPCParam->GetMWPCReadout()==kTRUE){
2056 fTPCParam->Transform2to2NearestWire(xyz,index);
2057 Float_t dx=xyz[0]-x1;
2058 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2060 //add nonisochronity (not implemented yet)
2064 ClassImp(AliTPCdigit)
2066 //_____________________________________________________________________________
2067 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2071 // Creates a TPC digit object
2073 fSector = digits[0];
2074 fPadRow = digits[1];
2077 fSignal = digits[4];
2083 //_____________________________________________________________________________
2084 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2088 // Creates a TPC hit object
2099 //________________________________________________________________________
2100 // Additional code because of the AliTPCTrackHits
2102 void AliTPC::MakeBranch2(Option_t *option,const char *file)
2105 // Create a new branch in the current Root Tree
2106 // The branch of fHits is automatically split
2107 // MI change 14.09.2000
2108 if (fHitType&2==0) return;
2109 char branchname[10];
2110 sprintf(branchname,"%s2",GetName());
2112 // Get the pointer to the header
2113 const char *cH = strstr(option,"H");
2115 if (fTrackHits && gAlice->TreeH() && cH) {
2116 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHits,
2117 gAlice->TreeH(),fBufferSize,1);
2118 gAlice->TreeH()->GetListOfBranches()->Add(branch);
2120 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2121 const char folder [] = "RunMC/Event/Data";
2123 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2124 Publish(folder,&fTrackHits,branchname);
2126 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2127 TDirectory *wd = gDirectory;
2129 TIter next( b->GetListOfBranches());
2130 while ((b=(TBranch*)next())) {
2135 cout << "Diverting branch " << branchname << " to file " << file << endl;
2140 void AliTPC::SetTreeAddress()
2142 if (fHitType&1) AliDetector::SetTreeAddress();
2143 if (fHitType&2) SetTreeAddress2();
2146 void AliTPC::SetTreeAddress2()
2149 // Set branch address for the TrackHits Tree
2152 char branchname[20];
2153 sprintf(branchname,"%s2",GetName());
2155 // Branch address for hit tree
2156 TTree *treeH = gAlice->TreeH();
2158 branch = treeH->GetBranch(branchname);
2159 if (branch) branch->SetAddress(&fTrackHits);
2163 void AliTPC::FinishPrimary()
2165 if (fTrackHits) fTrackHits->FlushHitStack();
2169 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2172 // add hit to the list
2175 int primary = gAlice->GetPrimary(track);
2176 gAlice->Particle(primary)->SetBit(kKeepBit);
2180 gAlice->FlagTrack(track);
2182 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2183 //if (hit->fTrack!=rtrack)
2184 // cout<<"bad track number\n";
2186 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2187 hits[1],hits[2],(Int_t)hits[3]);
2190 void AliTPC::ResetHits()
2192 if (fHitType&1) AliDetector::ResetHits();
2193 if (fHitType&2) ResetHits2();
2196 void AliTPC::ResetHits2()
2200 if (fTrackHits) fTrackHits->Clear();
2203 AliHit* AliTPC::FirstHit(Int_t track)
2205 if (fHitType&2) return FirstHit2(track);
2206 return AliDetector::FirstHit(track);
2208 AliHit* AliTPC::NextHit()
2210 if (fHitType&2) return NextHit2();
2211 return AliDetector::NextHit();
2214 AliHit* AliTPC::FirstHit2(Int_t track)
2217 // Initialise the hit iterator
2218 // Return the address of the first hit for track
2219 // If track>=0 the track is read from disk
2220 // while if track<0 the first hit of the current
2221 // track is returned
2224 gAlice->ResetHits();
2225 gAlice->TreeH()->GetEvent(track);
2229 fTrackHits->First();
2230 return fTrackHits->GetHit();
2235 AliHit* AliTPC::NextHit2()
2238 //Return the next hit for the current track
2242 return fTrackHits->GetHit();
2248 void AliTPC::LoadPoints(Int_t)
2252 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2255 if(fHitType==1) AliDetector::LoadPoints(a);
2256 else LoadPoints2(a);
2263 void AliTPC::RemapTrackHitIDs(Int_t *map)
2265 if (!fTrackHits) return;
2266 AliObjectArray * arr = fTrackHits->fTrackHitsInfo;
2267 for (UInt_t i=0;i<arr->GetSize();i++){
2268 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2269 info->fTrackID = map[info->fTrackID];
2275 //_____________________________________________________________________________
2276 void AliTPC::LoadPoints2(Int_t)
2279 // Store x, y, z of all hits in memory
2281 if (fTrackHits == 0) return;
2283 Int_t nhits = fTrackHits->GetEntriesFast();
2284 if (nhits == 0) return;
2285 Int_t tracks = gAlice->GetNtrack();
2286 if (fPoints == 0) fPoints = new TObjArray(tracks);
2289 Int_t *ntrk=new Int_t[tracks];
2290 Int_t *limi=new Int_t[tracks];
2291 Float_t **coor=new Float_t*[tracks];
2292 for(Int_t i=0;i<tracks;i++) {
2298 AliPoints *points = 0;
2301 Int_t chunk=nhits/4+1;
2303 // Loop over all the hits and store their position
2305 ahit = FirstHit2(-1);
2306 //for (Int_t hit=0;hit<nhits;hit++) {
2308 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2309 trk=ahit->GetTrack();
2310 if(ntrk[trk]==limi[trk]) {
2312 // Initialise a new track
2313 fp=new Float_t[3*(limi[trk]+chunk)];
2315 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2316 delete [] coor[trk];
2323 fp[3*ntrk[trk] ] = ahit->X();
2324 fp[3*ntrk[trk]+1] = ahit->Y();
2325 fp[3*ntrk[trk]+2] = ahit->Z();
2330 for(trk=0; trk<tracks; ++trk) {
2332 points = new AliPoints();
2333 points->SetMarkerColor(GetMarkerColor());
2334 points->SetMarkerSize(GetMarkerSize());
2335 points->SetDetector(this);
2336 points->SetParticle(trk);
2337 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2338 fPoints->AddAt(points,trk);
2339 delete [] coor[trk];
2349 //_____________________________________________________________________________
2350 void AliTPC::LoadPoints3(Int_t)
2353 // Store x, y, z of all hits in memory
2354 // - only intersection point with pad row
2355 if (fTrackHits == 0) return;
2357 Int_t nhits = fTrackHits->GetEntriesFast();
2358 if (nhits == 0) return;
2359 Int_t tracks = gAlice->GetNtrack();
2360 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2361 fPoints->Expand(2*tracks);
2364 Int_t *ntrk=new Int_t[tracks];
2365 Int_t *limi=new Int_t[tracks];
2366 Float_t **coor=new Float_t*[tracks];
2367 for(Int_t i=0;i<tracks;i++) {
2373 AliPoints *points = 0;
2376 Int_t chunk=nhits/4+1;
2378 // Loop over all the hits and store their position
2380 ahit = FirstHit2(-1);
2381 //for (Int_t hit=0;hit<nhits;hit++) {
2385 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2386 trk=ahit->GetTrack();
2387 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2388 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2389 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2390 if (currentrow!=lastrow){
2391 lastrow = currentrow;
2392 //later calculate intersection point
2393 if(ntrk[trk]==limi[trk]) {
2395 // Initialise a new track
2396 fp=new Float_t[3*(limi[trk]+chunk)];
2398 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2399 delete [] coor[trk];
2406 fp[3*ntrk[trk] ] = ahit->X();
2407 fp[3*ntrk[trk]+1] = ahit->Y();
2408 fp[3*ntrk[trk]+2] = ahit->Z();
2415 for(trk=0; trk<tracks; ++trk) {
2417 points = new AliPoints();
2418 points->SetMarkerColor(GetMarkerColor()+1);
2419 points->SetMarkerStyle(5);
2420 points->SetMarkerSize(0.2);
2421 points->SetDetector(this);
2422 points->SetParticle(trk);
2423 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2424 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2425 fPoints->AddAt(points,tracks+trk);
2426 delete [] coor[trk];
2437 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2441 //fill clones array with intersection of current point with the
2446 const Int_t kcmaxhits=30000;
2447 TVector * xxxx = new TVector(kcmaxhits*4);
2448 TVector & xxx = *xxxx;
2449 Int_t maxhits = kcmaxhits;
2452 AliTPChit * tpcHit=0;
2453 tpcHit = (AliTPChit*)FirstHit2(-1);
2454 Int_t currentIndex=0;
2455 Int_t lastrow=-1; //last writen row
2458 if (tpcHit==0) continue;
2459 sector=tpcHit->fSector; // sector number
2460 ipart=tpcHit->Track();
2464 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2465 Int_t index[3]={1,sector,0};
2466 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2467 if (currentrow<0) continue;
2468 if (lastrow<0) lastrow=currentrow;
2469 if (currentrow==lastrow){
2470 if ( currentIndex>=maxhits){
2472 xxx.ResizeTo(4*maxhits);
2474 xxx(currentIndex*4)=x[0];
2475 xxx(currentIndex*4+1)=x[1];
2476 xxx(currentIndex*4+2)=x[2];
2477 xxx(currentIndex*4+3)=tpcHit->fQ;
2481 if (currentIndex>2){
2493 for (Int_t index=0;index<currentIndex;index++){
2495 x=x2=x3=x4=xxx(index*4);
2503 sumy+=xxx(index*4+1);
2504 sumxy+=xxx(index*4+1)*x;
2505 sumx2y+=xxx(index*4+1)*x2;
2506 sumz+=xxx(index*4+2);
2507 sumxz+=xxx(index*4+2)*x;
2508 sumx2z+=xxx(index*4+2)*x2;
2509 sumq+=xxx(index*4+3);
2511 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2512 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2513 sumx2*(sumx*sumx3-sumx2*sumx2);
2515 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2516 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2517 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2518 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2520 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2521 sumx2*(sumx*sumx2y-sumx2*sumxy);
2522 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2523 sumx2*(sumx*sumx2z-sumx2*sumxz);
2525 Float_t y=detay/det+centralPad;
2526 Float_t z=detaz/det;
2527 Float_t by=detby/det; //y angle
2528 Float_t bz=detbz/det; //z angle
2529 sumy/=Float_t(currentIndex);
2530 sumz/=Float_t(currentIndex);
2532 AliComplexCluster cl;
2538 cl.fTracks[0]=ipart;
2540 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2541 if (row!=0) row->InsertCluster(&cl);
2544 } //end of calculating cluster for given row
2546 } // end of loop over hits
2550 //_______________________________________________________________________________
2551 void AliTPC::Digits2Reco(Int_t eventnumber)
2553 // empty for a time being