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.47 2001/11/19 10:25:34 kowal2
19 Nearest integer instead of integer when converting to ADC counts
21 Revision 1.46 2001/11/07 06:47:12 kowal2
24 Revision 1.45 2001/11/03 13:33:48 kowal2
25 Updated algorithms in Hits2SDigits, SDigits2Digits,
29 Revision 1.44 2001/08/30 09:28:48 hristov
30 TTree names are explicitly set via SetName(name) and then Write() is called
32 Revision 1.43 2001/07/28 12:02:54 hristov
33 Branch split level set to 99
35 Revision 1.42 2001/07/28 11:38:52 hristov
36 Loop variable declared once
38 Revision 1.41 2001/07/28 10:53:50 hristov
39 Digitisation done according to the general scheme (M.Ivanov)
41 Revision 1.40 2001/07/27 13:03:14 hristov
42 Default Branch split level set to 99
44 Revision 1.39 2001/07/26 09:09:34 kowal2
45 Hits2Reco method added
47 Revision 1.38 2001/07/20 14:32:43 kowal2
48 Processing of many events possible now
50 Revision 1.37 2001/06/12 07:17:18 kowal2
51 Hits2SDigits method implemented (summable digits)
53 Revision 1.36 2001/05/16 14:57:25 alibrary
54 New files for folders and Stack
56 Revision 1.35 2001/05/08 16:02:22 kowal2
57 Updated material specifications
59 Revision 1.34 2001/05/08 15:00:15 hristov
60 Corrections for tracking in arbitrary magnenetic field. Changes towards a concept of global Alice track. Back propagation of reconstructed tracks (Yu.Belikov)
62 Revision 1.33 2001/04/03 12:40:43 kowal2
65 Revision 1.32 2001/03/12 17:47:36 hristov
66 Changes needed on Sun with CC 5.0
68 Revision 1.31 2001/03/12 08:21:50 kowal2
69 Corrected C++ bug in the material definitions
71 Revision 1.30 2001/03/01 17:34:47 kowal2
72 Correction due to the accuracy problem
74 Revision 1.29 2001/02/28 16:34:40 kowal2
75 Protection against nonphysical values of the avalanche size,
78 Revision 1.28 2001/01/26 19:57:19 hristov
79 Major upgrade of AliRoot code
81 Revision 1.27 2001/01/13 17:29:33 kowal2
82 Sun compiler correction
84 Revision 1.26 2001/01/10 07:59:43 kowal2
85 Corrections to load points from the noncompressed hits.
87 Revision 1.25 2000/11/02 07:25:31 kowal2
88 Changes due to the new hit structure.
91 Revision 1.24 2000/10/05 16:06:09 kowal2
92 Forward declarations. Changes due to a new class AliComplexCluster.
94 Revision 1.23 2000/10/02 21:28:18 fca
95 Removal of useless dependecies via forward declarations
97 Revision 1.22 2000/07/10 20:57:39 hristov
98 Update of TPC code and macros by M.Kowalski
100 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
101 Changes to obey the coding rules
103 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
104 Splitted from AliTPCtracking
106 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
107 Changed parameter settings
109 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
111 Defaults loaded automatically (hard-wired)
112 Optional parameters can be set via macro called in the constructor
114 Revision 1.19 2000/04/18 19:00:59 fca
115 Small bug fixes to TPC files
117 Revision 1.18 2000/04/17 09:37:33 kowal2
118 removed obsolete AliTPCDigitsDisplay.C
120 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
122 New, experimental data structure from M. Ivanov
123 New tracking algorithm
124 Different pad geometry for different sectors
125 Digitization rewritten
127 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
128 Not used anymore - removed
130 Revision 1.17 2000/01/19 17:17:30 fca
131 Introducing a list of lists of hits -- more hits allowed for detector now
133 Revision 1.16 1999/11/05 09:29:23 fca
134 Accept only signals > 0
136 Revision 1.15 1999/10/08 06:26:53 fca
137 Removed ClustersIndex - not used anymore
139 Revision 1.14 1999/09/29 09:24:33 fca
140 Introduction of the Copyright and cvs Log
144 ///////////////////////////////////////////////////////////////////////////////
146 // Time Projection Chamber //
147 // This class contains the basic functions for the Time Projection Chamber //
148 // detector. Functions specific to one particular geometry are //
149 // contained in the derived classes //
153 <img src="picts/AliTPCClass.gif">
158 ///////////////////////////////////////////////////////////////////////////////
166 #include <TGeometry.h>
169 #include <TObjectTable.h>
170 #include "TParticle.h"
176 #include <iostream.h>
183 #include "AliTPCParamSR.h"
184 #include "AliTPCPRF2D.h"
185 #include "AliTPCRF1D.h"
186 #include "AliDigits.h"
187 #include "AliSimDigits.h"
188 #include "AliTPCTrackHits.h"
189 #include "AliPoints.h"
190 #include "AliArrayBranch.h"
193 #include "AliTPCDigitsArray.h"
194 #include "AliComplexCluster.h"
195 #include "AliClusters.h"
196 #include "AliTPCClustersRow.h"
197 #include "AliTPCClustersArray.h"
199 #include "AliTPCcluster.h"
200 #include "AliTPCclusterer.h"
201 #include "AliTPCtracker.h"
203 #include <TInterpreter.h>
210 //_____________________________________________________________________________
214 // Default constructor
231 //_____________________________________________________________________________
232 AliTPC::AliTPC(const char *name, const char *title)
233 : AliDetector(name,title)
236 // Standard constructor
240 // Initialise arrays of hits and digits
241 fHits = new TClonesArray("AliTPChit", 176);
242 gAlice->AddHitList(fHits);
248 fTrackHits = new AliTPCTrackHits; //MI - 13.09.2000
249 fTrackHits->SetHitPrecision(0.002);
250 fTrackHits->SetStepPrecision(0.003);
251 fTrackHits->SetMaxDistance(100);
256 // Initialise counters
262 // Initialise color attributes
263 SetMarkerColor(kYellow);
266 // Set TPC parameters
270 if (!strcmp(title,"Default")) {
271 fTPCParam = new AliTPCParamSR;
273 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
279 //_____________________________________________________________________________
290 delete fTrackHits; //MI 15.09.2000
291 if (fNoiseTable) delete [] fNoiseTable;
295 //_____________________________________________________________________________
296 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
299 // Add a hit to the list
301 // TClonesArray &lhits = *fHits;
302 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
304 TClonesArray &lhits = *fHits;
305 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
308 AddHit2(track,vol,hits);
311 //_____________________________________________________________________________
312 void AliTPC::BuildGeometry()
316 // Build TPC ROOT TNode geometry for the event display
321 const int kColorTPC=19;
322 char name[5], title[25];
323 const Double_t kDegrad=TMath::Pi()/180;
324 const Double_t kRaddeg=180./TMath::Pi();
327 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
328 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
330 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
331 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
333 Int_t nLo = fTPCParam->GetNInnerSector()/2;
334 Int_t nHi = fTPCParam->GetNOuterSector()/2;
336 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
337 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
338 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
339 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
342 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
343 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
349 // Get ALICE top node
352 nTop=gAlice->GetGeometry()->GetNode("alice");
356 rl = fTPCParam->GetInnerRadiusLow();
357 ru = fTPCParam->GetInnerRadiusUp();
361 sprintf(name,"LS%2.2d",i);
363 sprintf(title,"TPC low sector %3d",i);
366 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
367 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
368 tubs->SetNumberOfDivisions(1);
370 nNode = new TNode(name,title,name,0,0,0,"");
371 nNode->SetLineColor(kColorTPC);
377 rl = fTPCParam->GetOuterRadiusLow();
378 ru = fTPCParam->GetOuterRadiusUp();
381 sprintf(name,"US%2.2d",i);
383 sprintf(title,"TPC upper sector %d",i);
385 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
386 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
387 tubs->SetNumberOfDivisions(1);
389 nNode = new TNode(name,title,name,0,0,0,"");
390 nNode->SetLineColor(kColorTPC);
396 //_____________________________________________________________________________
397 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
400 // Calculate distance from TPC to mouse on the display
406 void AliTPC::Clusters2Tracks(TFile *of) {
407 //-----------------------------------------------------------------
408 // This is a track finder.
409 //-----------------------------------------------------------------
410 AliTPCtracker tracker(fTPCParam);
411 tracker.Clusters2Tracks(gFile,of);
414 //_____________________________________________________________________________
415 void AliTPC::CreateMaterials()
417 //-----------------------------------------------
418 // Create Materials for for TPC simulations
419 //-----------------------------------------------
421 //-----------------------------------------------------------------
422 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
423 //-----------------------------------------------------------------
425 Int_t iSXFLD=gAlice->Field()->Integ();
426 Float_t sXMGMX=gAlice->Field()->Max();
428 Float_t amat[5]; // atomic numbers
429 Float_t zmat[5]; // z
430 Float_t wmat[5]; // proportions
436 //***************** Gases *************************
438 //-------------------------------------------------
440 //-------------------------------------------------
451 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
461 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
464 //--------------------------------------------------------------
466 //--------------------------------------------------------------
483 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
485 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
500 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
502 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
518 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
520 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
522 //----------------------------------------------------------------
523 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
524 //----------------------------------------------------------------
530 Float_t rho,absl,X0,buf[1];
534 for(nc = 0;nc<fNoComp;nc++)
537 // retrive material constants
539 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
544 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
546 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
547 density += fMixtProp[nc]*rho; // density of the mixture
551 // mixture proportions by weight!
553 for(nc = 0;nc<fNoComp;nc++)
556 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
558 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
559 apure[nnc] : amol[nnc])/am;
563 // Drift gases 1 - nonsensitive, 2 - sensitive
565 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
566 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
575 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
578 //----------------------------------------------------------------------
580 //----------------------------------------------------------------------
602 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
624 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
642 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
660 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
662 // SiO2 - used later for the glass fiber
674 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
683 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
692 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
701 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
719 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
738 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
740 // Epoxy - C14 H20 O3
757 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
765 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
769 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,X0,absl,buf,nbuf);
773 wmat[0]=0.644; // by weight!
776 density=0.5*(1.25+2.265);
778 AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
782 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
784 wmat[0]=0.725; // by weight!
789 AliMixture(39,"G10",amat,zmat,density,2,wmat);
794 //----------------------------------------------------------
795 // tracking media for gases
796 //----------------------------------------------------------
798 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
799 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
800 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
801 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
803 //-----------------------------------------------------------
804 // tracking media for solids
805 //-----------------------------------------------------------
807 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
808 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
809 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
810 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
811 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
812 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
813 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
814 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
815 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
816 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
817 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
818 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
822 void AliTPC::GenerNoise(Int_t tablesize)
825 //Generate table with noise
832 if (fNoiseTable) delete[] fNoiseTable;
833 fNoiseTable = new Float_t[tablesize];
834 fNoiseDepth = tablesize;
835 fCurrentNoise =0; //!index of the noise in the noise table
837 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
838 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
841 Float_t AliTPC::GetNoise()
843 // get noise from table
844 // if ((fCurrentNoise%10)==0)
845 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
846 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
847 return fNoiseTable[fCurrentNoise++];
848 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
853 void AliTPC::Digits2Clusters(TFile *of, Int_t eventnumber)
855 //-----------------------------------------------------------------
856 // This is a simple cluster finder.
857 //-----------------------------------------------------------------
858 AliTPCclusterer::Digits2Clusters(fTPCParam,of,eventnumber);
861 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
862 extern Double_t SigmaZ2(Double_t, Double_t);
863 //_____________________________________________________________________________
864 void AliTPC::Hits2Clusters(TFile *of, Int_t eventn)
866 //--------------------------------------------------------
867 // TPC simple cluster generator from hits
868 // obtained from the TPC Fast Simulator
869 // The point errors are taken from the parametrization
870 //--------------------------------------------------------
872 //-----------------------------------------------------------------
873 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
874 //-----------------------------------------------------------------
875 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
876 // Jouri.Belikov@cern.ch
877 //----------------------------------------------------------------
879 /////////////////////////////////////////////////////////////////////////////
881 //---------------------------------------------------------------------
882 // ALICE TPC Cluster Parameters
883 //--------------------------------------------------------------------
887 // Cluster width in rphi
888 const Float_t kACrphi=0.18322;
889 const Float_t kBCrphi=0.59551e-3;
890 const Float_t kCCrphi=0.60952e-1;
891 // Cluster width in z
892 const Float_t kACz=0.19081;
893 const Float_t kBCz=0.55938e-3;
894 const Float_t kCCz=0.30428;
896 TDirectory *savedir=gDirectory;
899 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
903 //if(fDefaults == 0) SetDefaults();
905 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
907 TParticle *particle; // pointer to a given particle
908 AliTPChit *tpcHit; // pointer to a sigle TPC hit
912 Float_t pl,pt,tanth,rpad,ratio;
915 //---------------------------------------------------------------
916 // Get the access to the tracks
917 //---------------------------------------------------------------
919 TTree *tH = gAlice->TreeH();
920 Stat_t ntracks = tH->GetEntries();
922 //Switch to the output file
927 sprintf(cname,"TreeC_TPC_%d",eventn);
929 fTPCParam->Write(fTPCParam->GetTitle());
930 AliTPCClustersArray carray;
931 carray.Setup(fTPCParam);
932 carray.SetClusterType("AliTPCcluster");
935 Int_t nclusters=0; //cluster counter
937 //------------------------------------------------------------
938 // Loop over all sectors (72 sectors for 20 deg
939 // segmentation for both lower and upper sectors)
940 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
941 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
943 // First cluster for sector 0 starts at "0"
944 //------------------------------------------------------------
946 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
948 fTPCParam->AdjustCosSin(isec,cph,sph);
950 //------------------------------------------------------------
952 //------------------------------------------------------------
954 for(Int_t track=0;track<ntracks;track++){
958 // Get number of the TPC hits
960 // nhits=fHits->GetEntriesFast();
963 tpcHit = (AliTPChit*)FirstHit(-1);
967 // for(Int_t hit=0;hit<nhits;hit++){
968 //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
972 if (tpcHit->fQ == 0.) {
973 tpcHit = (AliTPChit*) NextHit();
974 continue; //information about track (I.Belikov)
976 sector=tpcHit->fSector; // sector number
979 // if(sector != isec) continue; //terminate iteration
982 tpcHit = (AliTPChit*) NextHit();
985 ipart=tpcHit->Track();
986 particle=gAlice->Particle(ipart);
989 if(pt < 1.e-9) pt=1.e-9;
991 tanth = TMath::Abs(tanth);
992 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
993 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
995 // space-point resolutions
997 sigmaRphi=SigmaY2(rpad,tanth,pt);
998 sigmaZ =SigmaZ2(rpad,tanth );
1002 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
1003 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
1005 // temporary protection
1007 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
1008 if(sigmaZ < 0.) sigmaZ=0.4e-3;
1009 if(clRphi < 0.) clRphi=2.5e-3;
1010 if(clZ < 0.) clZ=2.5e-5;
1015 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
1016 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1018 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
1019 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
1020 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
1021 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
1022 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1023 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
1024 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
1025 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
1026 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
1027 xyz[2]=sigmaRphi; // fSigmaY2
1028 xyz[3]=sigmaZ; // fSigmaZ2
1029 xyz[4]=tpcHit->fQ; // q
1031 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
1032 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
1034 Int_t tracks[3]={tpcHit->Track(), -1, -1};
1035 AliTPCcluster cluster(tracks,xyz);
1037 clrow->InsertCluster(&cluster); nclusters++;
1039 tpcHit = (AliTPChit*)NextHit();
1042 } // end of loop over hits
1044 } // end of loop over tracks
1046 Int_t nrows=fTPCParam->GetNRow(isec);
1047 for (Int_t irow=0; irow<nrows; irow++) {
1048 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
1049 if (!clrow) continue;
1050 carray.StoreRow(isec,irow);
1051 carray.ClearRow(isec,irow);
1054 } // end of loop over sectors
1056 cerr<<"Number of made clusters : "<<nclusters<<" \n";
1058 carray.GetTree()->SetName(cname);
1059 carray.GetTree()->Write();
1061 savedir->cd(); //switch back to the input file
1063 } // end of function
1065 //_________________________________________________________________
1066 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1068 //--------------------------------------------------------
1069 //calculate exact cross point of track and given pad row
1070 //resulting values are expressed in "digit" coordinata
1071 //--------------------------------------------------------
1073 //-----------------------------------------------------------------
1074 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1075 //-----------------------------------------------------------------
1077 if (fClustersArray==0){
1081 TParticle *particle; // pointer to a given particle
1082 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1083 // Int_t sector,nhits;
1085 const Int_t kcmaxhits=30000;
1086 TVector * xxxx = new TVector(kcmaxhits*4);
1087 TVector & xxx = *xxxx;
1088 Int_t maxhits = kcmaxhits;
1089 //construct array for each padrow
1090 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
1091 fClustersArray->CreateRow(isec,i);
1093 //---------------------------------------------------------------
1094 // Get the access to the tracks
1095 //---------------------------------------------------------------
1097 TTree *tH = gAlice->TreeH();
1098 Stat_t ntracks = tH->GetEntries();
1099 Int_t npart = gAlice->GetNtrack();
1100 TBranch * br = tH->GetBranch("fTrackHitsInfo");
1103 if (fHitType&2) branch = tH->GetBranch("TPC2");
1104 else branch = tH->GetBranch("TPC");
1106 //------------------------------------------------------------
1108 //------------------------------------------------------------
1110 for(Int_t track=0;track<ntracks;track++){
1111 Bool_t isInSector=kTRUE;
1116 br->GetEvent(track);
1117 AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
1118 for (UInt_t j=0;j<ar->GetSize();j++){
1119 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
1122 if (!isInSector) continue;
1124 branch->GetEntry(track); // get next track
1128 // Get number of the TPC hits and a pointer
1131 //nhits=fHits->GetEntriesFast();
1135 Int_t currentIndex=0;
1136 Int_t lastrow=-1; //last writen row
1140 tpcHit = (AliTPChit*)FirstHit(-1);
1143 Int_t sector=tpcHit->fSector; // sector number
1145 tpcHit = (AliTPChit*) NextHit();
1149 //for(Int_t hit=0;hit<nhits;hit++){
1150 //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
1151 //if (tpcHit==0) continue;
1152 //sector=tpcHit->fSector; // sector number
1153 //if(sector != isec) continue;
1154 ipart=tpcHit->Track();
1155 if (ipart<npart) particle=gAlice->Particle(ipart);
1159 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1160 Int_t index[3]={1,isec,0};
1161 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1162 if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
1163 if (lastrow<0) lastrow=currentrow;
1164 if (currentrow==lastrow){
1165 if ( currentIndex>=maxhits){
1167 xxx.ResizeTo(4*maxhits);
1169 xxx(currentIndex*4)=x[0];
1170 xxx(currentIndex*4+1)=x[1];
1171 xxx(currentIndex*4+2)=x[2];
1172 xxx(currentIndex*4+3)=tpcHit->fQ;
1176 if (currentIndex>2){
1188 for (Int_t index=0;index<currentIndex;index++){
1190 x=x2=x3=x4=xxx(index*4);
1198 sumy+=xxx(index*4+1);
1199 sumxy+=xxx(index*4+1)*x;
1200 sumx2y+=xxx(index*4+1)*x2;
1201 sumz+=xxx(index*4+2);
1202 sumxz+=xxx(index*4+2)*x;
1203 sumx2z+=xxx(index*4+2)*x2;
1204 sumq+=xxx(index*4+3);
1206 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1207 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1208 sumx2*(sumx*sumx3-sumx2*sumx2);
1210 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1211 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1212 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1213 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1215 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1216 sumx2*(sumx*sumx2y-sumx2*sumxy);
1217 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1218 sumx2*(sumx*sumx2z-sumx2*sumxz);
1220 Float_t y=detay/det+centralPad;
1221 Float_t z=detaz/det;
1222 Float_t by=detby/det; //y angle
1223 Float_t bz=detbz/det; //z angle
1224 sumy/=Float_t(currentIndex);
1225 sumz/=Float_t(currentIndex);
1227 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1229 AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
1230 // AliComplexCluster cl;
1236 cl->fTracks[0]=ipart;
1240 } //end of calculating cluster for given row
1243 tpcHit = (AliTPChit*)NextHit();
1244 } // end of loop over hits
1245 } // end of loop over tracks
1246 //write padrows to tree
1247 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1248 fClustersArray->StoreRow(isec,ii);
1249 fClustersArray->ClearRow(isec,ii);
1258 void AliTPC::SDigits2Digits2(Int_t eventnumber)
1260 //create digits from summable digits
1261 GenerNoise(500000); //create teble with noise
1264 sprintf(sname,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1265 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1267 //conect tree with sSDigits
1268 TTree *t = (TTree *)gDirectory->Get(sname);
1269 AliSimDigits digarr, *dummy=&digarr;
1270 t->GetBranch("Segment")->SetAddress(&dummy);
1271 Stat_t nentries = t->GetEntries();
1273 // set zero suppression
1275 fTPCParam->SetZeroSup(2);
1277 // get zero suppression
1279 Int_t zerosup = fTPCParam->GetZeroSup();
1281 //make tree with digits
1283 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1284 arr->SetClass("AliSimDigits");
1285 arr->Setup(fTPCParam);
1286 arr->MakeTree(fDigitsFile);
1288 AliTPCParam * par =fTPCParam;
1290 //Loop over segments of the TPC
1292 for (Int_t n=0; n<nentries; n++) {
1295 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1296 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
1299 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1300 Int_t nrows = digrow->GetNRows();
1301 Int_t ncols = digrow->GetNCols();
1303 digrow->ExpandBuffer();
1304 digarr.ExpandBuffer();
1305 digrow->ExpandTrackBuffer();
1306 digarr.ExpandTrackBuffer();
1309 Short_t * pamp0 = digarr.GetDigits();
1310 Int_t * ptracks0 = digarr.GetTracks();
1311 Short_t * pamp1 = digrow->GetDigits();
1312 Int_t * ptracks1 = digrow->GetTracks();
1313 Int_t nelems =nrows*ncols;
1314 Int_t saturation = fTPCParam->GetADCSat();
1315 //use internal structure of the AliDigits - for speed reason
1316 //if you cahnge implementation
1317 //of the Alidigits - it must be rewriten -
1318 for (Int_t i= 0; i<nelems; i++){
1319 // Float_t q = *pamp0;
1320 //q/=16.; //conversion faktor
1321 //Float_t noise= GetNoise();
1323 //q= TMath::Nint(q);
1324 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1326 if (q>saturation) q=saturation;
1328 //if (ptracks0[0]==0)
1331 ptracks1[0]=ptracks0[0];
1332 ptracks1[nelems]=ptracks0[nelems];
1333 ptracks1[2*nelems]=ptracks0[2*nelems];
1341 arr->StoreRow(sec,row);
1342 arr->ClearRow(sec,row);
1343 // cerr<<sec<<"\t"<<row<<"\n";
1350 arr->GetTree()->SetName(dname);
1351 arr->GetTree()->Write();
1354 //_________________________________________
1355 void AliTPC::Merge(TTree * intree, Int_t *mask, Int_t nin, Int_t outid)
1358 //intree - pointer to array of trees with s digits
1359 //mask - mask for each
1360 //nin - number of inputs
1361 //outid - event number of the output event
1363 //create digits from summable digits
1364 //conect tree with sSDigits
1367 AliSimDigits ** digarr = new AliSimDigits*[nin];
1368 for (Int_t i1=0;i1<nin; i1++){
1370 intree[i1].GetBranch("Segment")->SetAddress(&digarr[i1]);
1372 Stat_t nentries = intree[0].GetEntries();
1374 //make tree with digits
1376 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),outid);
1377 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1378 arr->SetClass("AliSimDigits");
1379 arr->Setup(fTPCParam);
1380 arr->MakeTree(fDigitsFile);
1382 // set zero suppression
1384 fTPCParam->SetZeroSup(2);
1386 // get zero suppression
1388 Int_t zerosup = fTPCParam->GetZeroSup();
1391 AliTPCParam * par =fTPCParam;
1393 //Loop over segments of the TPC
1394 for (Int_t n=0; n<nentries; n++) {
1396 for (Int_t i=0;i<nin; i++){
1397 //connect proper digits
1398 intree[i].GetEvent(n);
1399 digarr[i]->ExpandBuffer();
1400 digarr[i]->ExpandTrackBuffer();
1403 if (!par->AdjustSectorRow(digarr[0]->GetID(),sec,row)) {
1404 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr[0]->GetID()<<endl;
1408 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1409 Int_t nrows = digrow->GetNRows();
1410 Int_t ncols = digrow->GetNCols();
1412 digrow->ExpandBuffer();
1413 digrow->ExpandTrackBuffer();
1415 for (Int_t rows=0;rows<nrows; rows++){
1416 for (Int_t col=0;col<ncols; col++){
1418 Int_t label[1000]; // i hope no more than 300 events merged
1420 // looop over digits
1421 for (Int_t i=0;i<nin; i++){
1422 q += digarr[i]->GetDigitFast(rows,col);
1423 for (Int_t tr=0;tr<3;tr++) {
1424 Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
1426 label[labptr]=lab+mask[i]-2;
1432 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac()*16);
1434 q/=16; //conversion faktor
1439 if(q > fTPCParam->GetADCSat()) q = (Short_t)(fTPCParam->GetADCSat());
1440 digrow->SetDigitFast(q,rows,col);
1441 for (Int_t tr=0;tr<3;tr++){
1443 ((AliSimDigits*)digrow)->SetTrackIDFast(label[tr],rows,col,tr);
1445 ((AliSimDigits*)digrow)->SetTrackIDFast(-1,rows,col,tr);
1450 arr->StoreRow(sec,row);
1452 arr->ClearRow(sec,row);
1453 // cerr<<sec<<"\t"<<row<<"\n";
1457 arr->GetTree()->SetName(dname);
1458 arr->GetTree()->Write();
1464 /*_________________________________________
1465 void AliTPC::SDigits2Digits(Int_t eventnumber)
1469 cerr<<"Digitizing TPC...\n";
1471 Hits2Digits(eventnumber);
1476 // char treeName[100];
1478 // sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1480 // GetDigitsArray()->GetTree()->Write(treeName);
1483 //__________________________________________________________________
1484 void AliTPC::SetDefaults(){
1487 cerr<<"Setting default parameters...\n";
1489 // Set response functions
1491 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1492 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1493 AliTPCPRF2D * prfouter = new AliTPCPRF2D;
1494 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1495 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1496 rf->SetOffset(3*param->GetZSigma());
1499 TDirectory *savedir=gDirectory;
1500 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1502 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1505 prfinner->Read("prf_07504_Gati_056068_d02");
1506 prfouter->Read("prf_10006_Gati_047051_d03");
1510 param->SetInnerPRF(prfinner);
1511 param->SetOuterPRF(prfouter);
1512 param->SetTimeRF(rf);
1522 //__________________________________________________________________
1523 void AliTPC::Hits2Digits(Int_t eventnumber)
1525 //----------------------------------------------------
1526 // Loop over all sectors for a single event
1527 //----------------------------------------------------
1530 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1531 GenerNoise(500000); //create teble with noise
1533 //setup TPCDigitsArray
1535 if(GetDigitsArray()) delete GetDigitsArray();
1537 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1538 arr->SetClass("AliSimDigits");
1539 arr->Setup(fTPCParam);
1540 arr->MakeTree(fDigitsFile);
1541 SetDigitsArray(arr);
1543 fDigitsSwitch=0; // standard digits
1545 cerr<<"Digitizing TPC...\n";
1547 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1553 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1555 GetDigitsArray()->GetTree()->SetName(treeName);
1556 GetDigitsArray()->GetTree()->Write();
1563 //__________________________________________________________________
1564 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1567 //-----------------------------------------------------------
1568 // summable digits - 16 bit "ADC", no noise, no saturation
1569 //-----------------------------------------------------------
1571 //----------------------------------------------------
1572 // Loop over all sectors for a single event
1573 //----------------------------------------------------
1576 if(fDefaults == 0) SetDefaults();
1577 GenerNoise(500000); //create table with noise
1578 //setup TPCDigitsArray
1580 if(GetDigitsArray()) delete GetDigitsArray();
1582 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1583 arr->SetClass("AliSimDigits");
1584 arr->Setup(fTPCParam);
1585 arr->MakeTree(fDigitsFile);
1586 SetDigitsArray(arr);
1588 cerr<<"Digitizing TPC...\n";
1590 fDigitsSwitch=1; // summable digits
1592 // set zero suppression to "0"
1594 fTPCParam->SetZeroSup(0);
1596 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1603 sprintf(treeName,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1605 GetDigitsArray()->GetTree()->SetName(treeName);
1606 GetDigitsArray()->GetTree()->Write();
1612 //__________________________________________________________________
1613 void AliTPC::Hits2SDigits()
1616 //-----------------------------------------------------------
1617 // summable digits - 16 bit "ADC", no noise, no saturation
1618 //-----------------------------------------------------------
1620 //----------------------------------------------------
1621 // Loop over all sectors for a single event
1622 //----------------------------------------------------
1623 //MI change - for pp run
1624 Int_t eventnumber = gAlice->GetEvNumber();
1626 if(fDefaults == 0) SetDefaults();
1628 //setup TPCDigitsArray
1630 if(GetDigitsArray()) delete GetDigitsArray();
1632 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1633 arr->SetClass("AliSimDigits");
1634 arr->Setup(fTPCParam);
1635 arr->MakeTree(fDigitsFile);
1636 SetDigitsArray(arr);
1638 cerr<<"Digitizing TPC...\n";
1640 // fDigitsSwitch=1; // summable digits -for the moment direct
1642 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1648 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1650 GetDigitsArray()->GetTree()->SetName(treeName);
1651 GetDigitsArray()->GetTree()->Write();
1656 //_____________________________________________________________________________
1657 void AliTPC::Hits2DigitsSector(Int_t isec)
1659 //-------------------------------------------------------------------
1660 // TPC conversion from hits to digits.
1661 //-------------------------------------------------------------------
1663 //-----------------------------------------------------------------
1664 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1665 //-----------------------------------------------------------------
1667 //-------------------------------------------------------
1668 // Get the access to the track hits
1669 //-------------------------------------------------------
1671 // check if the parameters are set - important if one calls this method
1672 // directly, not from the Hits2Digits
1674 if(fDefaults == 0) SetDefaults();
1676 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1677 Stat_t ntracks = tH->GetEntries();
1681 //-------------------------------------------
1682 // Only if there are any tracks...
1683 //-------------------------------------------
1687 //printf("*** Processing sector number %d ***\n",isec);
1689 Int_t nrows =fTPCParam->GetNRow(isec);
1691 row= new TObjArray* [nrows];
1693 MakeSector(isec,nrows,tH,ntracks,row);
1695 //--------------------------------------------------------
1696 // Digitize this sector, row by row
1697 // row[i] is the pointer to the TObjArray of TVectors,
1698 // each one containing electrons accepted on this
1699 // row, assigned into tracks
1700 //--------------------------------------------------------
1704 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1706 for (i=0;i<nrows;i++){
1708 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1710 DigitizeRow(i,isec,row);
1712 fDigitsArray->StoreRow(isec,i);
1714 Int_t ndig = dig->GetDigitSize();
1717 //printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1719 fDigitsArray->ClearRow(isec,i);
1722 } // end of the sector digitization
1724 for(i=0;i<nrows;i++){
1729 delete [] row; // delete the array of pointers to TObjArray-s
1733 } // end of Hits2DigitsSector
1736 //_____________________________________________________________________________
1737 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1739 //-----------------------------------------------------------
1740 // Single row digitization, coupling from the neighbouring
1741 // rows taken into account
1742 //-----------------------------------------------------------
1744 //-----------------------------------------------------------------
1745 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1746 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1747 //-----------------------------------------------------------------
1750 Float_t zerosup = fTPCParam->GetZeroSup();
1751 Int_t nrows =fTPCParam->GetNRow(isec);
1752 fCurrentIndex[1]= isec;
1755 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1756 Int_t nofTbins = fTPCParam->GetMaxTBin();
1757 Int_t indexRange[4];
1759 // Integrated signal for this row
1760 // and a single track signal
1762 TMatrix *m1 = new TMatrix(0,nofPads,0,nofTbins); // integrated
1763 TMatrix *m2 = new TMatrix(0,nofPads,0,nofTbins); // single
1765 TMatrix &total = *m1;
1767 // Array of pointers to the label-signal list
1769 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1770 Float_t **pList = new Float_t* [nofDigits];
1774 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1778 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1779 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1780 for (Int_t row= row1;row<=row2;row++){
1781 Int_t nTracks= rows[row]->GetEntries();
1782 for (i1=0;i1<nTracks;i1++){
1783 fCurrentIndex[2]= row;
1784 fCurrentIndex[3]=irow;
1786 m2->Zero(); // clear single track signal matrix
1787 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1788 GetList(trackLabel,nofPads,m2,indexRange,pList);
1790 else GetSignal(rows[row],i1,0,m1,indexRange);
1796 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1797 for(Int_t ip=0;ip<nofPads;ip++){
1798 for(Int_t it=0;it<nofTbins;it++){
1800 Float_t q = total(ip,it);
1802 Int_t gi =it*nofPads+ip; // global index
1804 if(fDigitsSwitch == 0){
1806 // q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
1810 if(q <=zerosup) continue; // do not fill zeros
1811 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1816 if(q <= 0.) continue; // do not fill zeros
1817 if(q>2000.) q=2000.;
1823 // "real" signal or electronic noise (list = -1)?
1826 for(Int_t j1=0;j1<3;j1++){
1827 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1832 <A NAME="AliDigits"></A>
1833 using of AliDigits object
1836 dig->SetDigitFast((Short_t)q,it,ip);
1837 if (fDigitsArray->IsSimulated())
1839 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1840 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1841 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1845 } // end of loop over time buckets
1846 } // end of lop over pads
1849 // This row has been digitized, delete nonused stuff
1852 for(lp=0;lp<nofDigits;lp++){
1853 if(pList[lp]) delete [] pList[lp];
1862 } // end of DigitizeRow
1864 //_____________________________________________________________________________
1866 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, TMatrix *m1, TMatrix *m2,
1870 //---------------------------------------------------------------
1871 // Calculates 2-D signal (pad,time) for a single track,
1872 // returns a pointer to the signal matrix and the track label
1873 // No digitization is performed at this level!!!
1874 //---------------------------------------------------------------
1876 //-----------------------------------------------------------------
1877 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1878 // Modified: Marian Ivanov
1879 //-----------------------------------------------------------------
1883 tv = (TVector*)p1->At(ntr); // pointer to a track
1886 Float_t label = v(0);
1887 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1889 Int_t nElectrons = (tv->GetNrows()-1)/4;
1890 indexRange[0]=9999; // min pad
1891 indexRange[1]=-1; // max pad
1892 indexRange[2]=9999; //min time
1893 indexRange[3]=-1; // max time
1895 // Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1897 TMatrix &signal = *m1;
1898 TMatrix &total = *m2;
1900 // Loop over all electrons
1902 for(Int_t nel=0; nel<nElectrons; nel++){
1904 Float_t aval = v(idx+4);
1905 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1906 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1907 Int_t n = fTPCParam->CalcResponse(xyz,fCurrentIndex,fCurrentIndex[3]);
1909 if (n>0) for (Int_t i =0; i<n; i++){
1910 Int_t *index = fTPCParam->GetResBin(i);
1911 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1912 if ( ( pad<(fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]))) && (pad>0)) {
1913 Int_t time=index[2];
1914 Float_t weight = fTPCParam->GetResWeight(i); //we normalise response to ADC channel
1915 weight *= eltoadcfac;
1917 if (m1!=0) signal(pad,time)+=weight;
1918 total(pad,time)+=weight;
1919 indexRange[0]=TMath::Min(indexRange[0],pad);
1920 indexRange[1]=TMath::Max(indexRange[1],pad);
1921 indexRange[2]=TMath::Min(indexRange[2],time);
1922 indexRange[3]=TMath::Max(indexRange[3],time);
1925 } // end of loop over electrons
1927 return label; // returns track label when finished
1930 //_____________________________________________________________________________
1931 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *indexRange,
1934 //----------------------------------------------------------------------
1935 // Updates the list of tracks contributing to digits for a given row
1936 //----------------------------------------------------------------------
1938 //-----------------------------------------------------------------
1939 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1940 //-----------------------------------------------------------------
1942 TMatrix &signal = *m;
1944 // lop over nonzero digits
1946 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1947 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1950 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1952 if(signal(ip,it)<0.5) continue;
1955 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1957 if(!pList[globalIndex]){
1960 // Create new list (6 elements - 3 signals and 3 labels),
1963 pList[globalIndex] = new Float_t [6];
1967 *pList[globalIndex] = -1.;
1968 *(pList[globalIndex]+1) = -1.;
1969 *(pList[globalIndex]+2) = -1.;
1970 *(pList[globalIndex]+3) = -1.;
1971 *(pList[globalIndex]+4) = -1.;
1972 *(pList[globalIndex]+5) = -1.;
1975 *pList[globalIndex] = label;
1976 *(pList[globalIndex]+3) = signal(ip,it);
1980 // check the signal magnitude
1982 Float_t highest = *(pList[globalIndex]+3);
1983 Float_t middle = *(pList[globalIndex]+4);
1984 Float_t lowest = *(pList[globalIndex]+5);
1987 // compare the new signal with already existing list
1990 if(signal(ip,it)<lowest) continue; // neglect this track
1994 if (signal(ip,it)>highest){
1995 *(pList[globalIndex]+5) = middle;
1996 *(pList[globalIndex]+4) = highest;
1997 *(pList[globalIndex]+3) = signal(ip,it);
1999 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2000 *(pList[globalIndex]+1) = *pList[globalIndex];
2001 *pList[globalIndex] = label;
2003 else if (signal(ip,it)>middle){
2004 *(pList[globalIndex]+5) = middle;
2005 *(pList[globalIndex]+4) = signal(ip,it);
2007 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2008 *(pList[globalIndex]+1) = label;
2011 *(pList[globalIndex]+5) = signal(ip,it);
2012 *(pList[globalIndex]+2) = label;
2016 } // end of loop over pads
2017 } // end of loop over time bins
2022 //___________________________________________________________________
2023 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
2024 Stat_t ntracks,TObjArray **row)
2027 //-----------------------------------------------------------------
2028 // Prepares the sector digitization, creates the vectors of
2029 // tracks for each row of this sector. The track vector
2030 // contains the track label and the position of electrons.
2031 //-----------------------------------------------------------------
2033 //-----------------------------------------------------------------
2034 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2035 //-----------------------------------------------------------------
2037 Float_t gasgain = fTPCParam->GetGasGain();
2041 AliTPChit *tpcHit; // pointer to a sigle TPC hit
2044 if (fHitType&2) branch = TH->GetBranch("TPC2");
2045 else branch = TH->GetBranch("TPC");
2048 //----------------------------------------------
2049 // Create TObjArray-s, one for each row,
2050 // each TObjArray will store the TVectors
2051 // of electrons, one TVector per each track.
2052 //----------------------------------------------
2054 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
2055 TVector **tracks = new TVector* [nrows]; //pointers to the track vectors
2056 for(i=0; i<nrows; i++){
2057 row[i] = new TObjArray;
2064 //--------------------------------------------------------------------
2065 // Loop over tracks, the "track" contains the full history
2066 //--------------------------------------------------------------------
2068 Int_t previousTrack,currentTrack;
2069 previousTrack = -1; // nothing to store so far!
2071 for(Int_t track=0;track<ntracks;track++){
2072 Bool_t isInSector=kTRUE;
2077 TBranch * br = TH->GetBranch("fTrackHitsInfo");
2078 br->GetEvent(track);
2079 AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
2080 for (UInt_t j=0;j<ar->GetSize();j++){
2081 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
2084 if (!isInSector) continue;
2086 branch->GetEntry(track); // get next track
2090 tpcHit = (AliTPChit*)FirstHit(-1);
2092 //--------------------------------------------------------------
2094 //--------------------------------------------------------------
2099 Int_t sector=tpcHit->fSector; // sector number
2100 // if(sector != isec) continue;
2102 tpcHit = (AliTPChit*) NextHit();
2106 currentTrack = tpcHit->Track(); // track number
2109 if(currentTrack != previousTrack){
2111 // store already filled fTrack
2113 for(i=0;i<nrows;i++){
2114 if(previousTrack != -1){
2115 if(nofElectrons[i]>0){
2116 TVector &v = *tracks[i];
2117 v(0) = previousTrack;
2118 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2119 row[i]->Add(tracks[i]);
2122 delete tracks[i]; // delete empty TVector
2128 tracks[i] = new TVector(481); // TVectors for the next fTrack
2130 } // end of loop over rows
2132 previousTrack=currentTrack; // update track label
2135 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2137 //---------------------------------------------------
2138 // Calculate the electron attachment probability
2139 //---------------------------------------------------
2142 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
2143 /fTPCParam->GetDriftV();
2145 Float_t attProb = fTPCParam->GetAttCoef()*
2146 fTPCParam->GetOxyCont()*time; // fraction!
2148 //-----------------------------------------------
2149 // Loop over electrons
2150 //-----------------------------------------------
2153 for(Int_t nel=0;nel<qI;nel++){
2154 // skip if electron lost due to the attachment
2155 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2160 // protection for the nonphysical avalanche size (10**6 maximum)
2162 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2163 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2166 TransportElectron(xyz,index); //MI change -august
2168 fTPCParam->GetPadRow(xyz,index); //MI change august
2169 rowNumber = index[2];
2170 //transform position to local digit coordinates
2171 //relative to nearest pad row
2172 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
2173 nofElectrons[rowNumber]++;
2174 //----------------------------------
2175 // Expand vector if necessary
2176 //----------------------------------
2177 if(nofElectrons[rowNumber]>120){
2178 Int_t range = tracks[rowNumber]->GetNrows();
2179 if((nofElectrons[rowNumber])>(range-1)/4){
2181 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
2185 TVector &v = *tracks[rowNumber];
2186 Int_t idx = 4*nofElectrons[rowNumber]-3;
2188 v(idx)= xyz[0]; // X - pad row coordinate
2189 v(idx+1)=xyz[1]; // Y - pad coordinate (along the pad-row)
2190 v(idx+2)=xyz[2]; // Z - time bin coordinate
2191 v(idx+3)=xyz[3]; // avalanche size
2192 } // end of loop over electrons
2194 tpcHit = (AliTPChit*)NextHit();
2196 } // end of loop over hits
2197 } // end of loop over tracks
2200 // store remaining track (the last one) if not empty
2203 for(i=0;i<nrows;i++){
2204 if(nofElectrons[i]>0){
2205 TVector &v = *tracks[i];
2206 v(0) = previousTrack;
2207 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2208 row[i]->Add(tracks[i]);
2217 delete [] nofElectrons;
2220 } // end of MakeSector
2223 //_____________________________________________________________________________
2227 // Initialise TPC detector after definition of geometry
2232 printf("\n%s: ",ClassName());
2233 for(i=0;i<35;i++) printf("*");
2234 printf(" TPC_INIT ");
2235 for(i=0;i<35;i++) printf("*");
2236 printf("\n%s: ",ClassName());
2238 for(i=0;i<80;i++) printf("*");
2243 //_____________________________________________________________________________
2244 void AliTPC::MakeBranch(Option_t* option, const char *file)
2247 // Create Tree branches for the TPC.
2249 Int_t buffersize = 4000;
2250 char branchname[10];
2251 sprintf(branchname,"%s",GetName());
2253 AliDetector::MakeBranch(option,file);
2255 const char *d = strstr(option,"D");
2257 if (fDigits && gAlice->TreeD() && d) {
2258 MakeBranchInTree(gAlice->TreeD(),
2259 branchname, &fDigits, buffersize, file);
2262 if (fHitType&2) MakeBranch2(option,file); // MI change 14.09.2000
2265 //_____________________________________________________________________________
2266 void AliTPC::ResetDigits()
2269 // Reset number of digits and the digits array for this detector
2272 if (fDigits) fDigits->Clear();
2275 //_____________________________________________________________________________
2276 void AliTPC::SetSecAL(Int_t sec)
2278 //---------------------------------------------------
2279 // Activate/deactivate selection for lower sectors
2280 //---------------------------------------------------
2282 //-----------------------------------------------------------------
2283 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2284 //-----------------------------------------------------------------
2289 //_____________________________________________________________________________
2290 void AliTPC::SetSecAU(Int_t sec)
2292 //----------------------------------------------------
2293 // Activate/deactivate selection for upper sectors
2294 //---------------------------------------------------
2296 //-----------------------------------------------------------------
2297 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2298 //-----------------------------------------------------------------
2303 //_____________________________________________________________________________
2304 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2306 //----------------------------------------
2307 // Select active lower sectors
2308 //----------------------------------------
2310 //-----------------------------------------------------------------
2311 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2312 //-----------------------------------------------------------------
2322 //_____________________________________________________________________________
2323 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2324 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2325 Int_t s11 , Int_t s12)
2327 //--------------------------------
2328 // Select active upper sectors
2329 //--------------------------------
2331 //-----------------------------------------------------------------
2332 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2333 //-----------------------------------------------------------------
2349 //_____________________________________________________________________________
2350 void AliTPC::SetSens(Int_t sens)
2353 //-------------------------------------------------------------
2354 // Activates/deactivates the sensitive strips at the center of
2355 // the pad row -- this is for the space-point resolution calculations
2356 //-------------------------------------------------------------
2358 //-----------------------------------------------------------------
2359 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2360 //-----------------------------------------------------------------
2366 void AliTPC::SetSide(Float_t side=0.)
2368 // choice of the TPC side
2373 //____________________________________________________________________________
2374 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2375 Float_t p2,Float_t p3)
2378 // gax mixture definition
2392 //_____________________________________________________________________________
2394 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2397 // electron transport taking into account:
2399 // 2.ExB at the wires
2400 // 3. nonisochronity
2402 // xyz and index must be already transformed to system 1
2405 fTPCParam->Transform1to2(xyz,index);
2408 Float_t driftl=xyz[2];
2409 if(driftl<0.01) driftl=0.01;
2410 driftl=TMath::Sqrt(driftl);
2411 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2412 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2413 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2414 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2415 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2419 if (fTPCParam->GetMWPCReadout()==kTRUE){
2421 fTPCParam->Transform2to2NearestWire(xyz,index);
2422 Float_t dx=xyz[0]-x1;
2423 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2425 //add nonisochronity (not implemented yet)
2429 ClassImp(AliTPCdigit)
2431 //_____________________________________________________________________________
2432 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2436 // Creates a TPC digit object
2438 fSector = digits[0];
2439 fPadRow = digits[1];
2442 fSignal = digits[4];
2448 //_____________________________________________________________________________
2449 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2453 // Creates a TPC hit object
2464 //________________________________________________________________________
2465 // Additional code because of the AliTPCTrackHits
2467 void AliTPC::MakeBranch2(Option_t *option,const char *file)
2470 // Create a new branch in the current Root Tree
2471 // The branch of fHits is automatically split
2472 // MI change 14.09.2000
2473 if (fHitType&2==0) return;
2474 char branchname[10];
2475 sprintf(branchname,"%s2",GetName());
2477 // Get the pointer to the header
2478 const char *cH = strstr(option,"H");
2480 if (fTrackHits && gAlice->TreeH() && cH) {
2481 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHits,
2482 gAlice->TreeH(),fBufferSize,99);
2483 gAlice->TreeH()->GetListOfBranches()->Add(branch);
2485 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2486 const char folder [] = "RunMC/Event/Data";
2488 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2489 Publish(folder,&fTrackHits,branchname);
2491 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2492 TDirectory *wd = gDirectory;
2494 TIter next( b->GetListOfBranches());
2495 while ((b=(TBranch*)next())) {
2500 cout << "Diverting branch " << branchname << " to file " << file << endl;
2505 void AliTPC::SetTreeAddress()
2507 if (fHitType&1) AliDetector::SetTreeAddress();
2508 if (fHitType&2) SetTreeAddress2();
2511 void AliTPC::SetTreeAddress2()
2514 // Set branch address for the TrackHits Tree
2517 char branchname[20];
2518 sprintf(branchname,"%s2",GetName());
2520 // Branch address for hit tree
2521 TTree *treeH = gAlice->TreeH();
2523 branch = treeH->GetBranch(branchname);
2524 if (branch) branch->SetAddress(&fTrackHits);
2528 void AliTPC::FinishPrimary()
2530 if (fTrackHits) fTrackHits->FlushHitStack();
2534 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2537 // add hit to the list
2540 int primary = gAlice->GetPrimary(track);
2541 gAlice->Particle(primary)->SetBit(kKeepBit);
2545 gAlice->FlagTrack(track);
2547 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2548 //if (hit->fTrack!=rtrack)
2549 // cout<<"bad track number\n";
2551 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2552 hits[1],hits[2],(Int_t)hits[3]);
2555 void AliTPC::ResetHits()
2557 if (fHitType&1) AliDetector::ResetHits();
2558 if (fHitType&2) ResetHits2();
2561 void AliTPC::ResetHits2()
2565 if (fTrackHits) fTrackHits->Clear();
2568 AliHit* AliTPC::FirstHit(Int_t track)
2570 if (fHitType&2) return FirstHit2(track);
2571 return AliDetector::FirstHit(track);
2573 AliHit* AliTPC::NextHit()
2575 if (fHitType&2) return NextHit2();
2576 return AliDetector::NextHit();
2579 AliHit* AliTPC::FirstHit2(Int_t track)
2582 // Initialise the hit iterator
2583 // Return the address of the first hit for track
2584 // If track>=0 the track is read from disk
2585 // while if track<0 the first hit of the current
2586 // track is returned
2589 gAlice->ResetHits();
2590 gAlice->TreeH()->GetEvent(track);
2594 fTrackHits->First();
2595 return fTrackHits->GetHit();
2600 AliHit* AliTPC::NextHit2()
2603 //Return the next hit for the current track
2607 return fTrackHits->GetHit();
2613 void AliTPC::LoadPoints(Int_t)
2617 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2620 if(fHitType==1) AliDetector::LoadPoints(a);
2621 else LoadPoints2(a);
2628 void AliTPC::RemapTrackHitIDs(Int_t *map)
2630 if (!fTrackHits) return;
2631 AliObjectArray * arr = fTrackHits->fTrackHitsInfo;
2632 for (UInt_t i=0;i<arr->GetSize();i++){
2633 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2634 info->fTrackID = map[info->fTrackID];
2640 //_____________________________________________________________________________
2641 void AliTPC::LoadPoints2(Int_t)
2644 // Store x, y, z of all hits in memory
2646 if (fTrackHits == 0) return;
2648 Int_t nhits = fTrackHits->GetEntriesFast();
2649 if (nhits == 0) return;
2650 Int_t tracks = gAlice->GetNtrack();
2651 if (fPoints == 0) fPoints = new TObjArray(tracks);
2654 Int_t *ntrk=new Int_t[tracks];
2655 Int_t *limi=new Int_t[tracks];
2656 Float_t **coor=new Float_t*[tracks];
2657 for(Int_t i=0;i<tracks;i++) {
2663 AliPoints *points = 0;
2666 Int_t chunk=nhits/4+1;
2668 // Loop over all the hits and store their position
2670 ahit = FirstHit2(-1);
2671 //for (Int_t hit=0;hit<nhits;hit++) {
2673 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2674 trk=ahit->GetTrack();
2675 if(ntrk[trk]==limi[trk]) {
2677 // Initialise a new track
2678 fp=new Float_t[3*(limi[trk]+chunk)];
2680 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2681 delete [] coor[trk];
2688 fp[3*ntrk[trk] ] = ahit->X();
2689 fp[3*ntrk[trk]+1] = ahit->Y();
2690 fp[3*ntrk[trk]+2] = ahit->Z();
2695 for(trk=0; trk<tracks; ++trk) {
2697 points = new AliPoints();
2698 points->SetMarkerColor(GetMarkerColor());
2699 points->SetMarkerSize(GetMarkerSize());
2700 points->SetDetector(this);
2701 points->SetParticle(trk);
2702 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2703 fPoints->AddAt(points,trk);
2704 delete [] coor[trk];
2714 //_____________________________________________________________________________
2715 void AliTPC::LoadPoints3(Int_t)
2718 // Store x, y, z of all hits in memory
2719 // - only intersection point with pad row
2720 if (fTrackHits == 0) return;
2722 Int_t nhits = fTrackHits->GetEntriesFast();
2723 if (nhits == 0) return;
2724 Int_t tracks = gAlice->GetNtrack();
2725 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2726 fPoints->Expand(2*tracks);
2729 Int_t *ntrk=new Int_t[tracks];
2730 Int_t *limi=new Int_t[tracks];
2731 Float_t **coor=new Float_t*[tracks];
2732 for(Int_t i=0;i<tracks;i++) {
2738 AliPoints *points = 0;
2741 Int_t chunk=nhits/4+1;
2743 // Loop over all the hits and store their position
2745 ahit = FirstHit2(-1);
2746 //for (Int_t hit=0;hit<nhits;hit++) {
2750 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2751 trk=ahit->GetTrack();
2752 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2753 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2754 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2755 if (currentrow!=lastrow){
2756 lastrow = currentrow;
2757 //later calculate intersection point
2758 if(ntrk[trk]==limi[trk]) {
2760 // Initialise a new track
2761 fp=new Float_t[3*(limi[trk]+chunk)];
2763 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2764 delete [] coor[trk];
2771 fp[3*ntrk[trk] ] = ahit->X();
2772 fp[3*ntrk[trk]+1] = ahit->Y();
2773 fp[3*ntrk[trk]+2] = ahit->Z();
2780 for(trk=0; trk<tracks; ++trk) {
2782 points = new AliPoints();
2783 points->SetMarkerColor(GetMarkerColor()+1);
2784 points->SetMarkerStyle(5);
2785 points->SetMarkerSize(0.2);
2786 points->SetDetector(this);
2787 points->SetParticle(trk);
2788 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2789 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2790 fPoints->AddAt(points,tracks+trk);
2791 delete [] coor[trk];
2802 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2806 //fill clones array with intersection of current point with the
2811 const Int_t kcmaxhits=30000;
2812 TVector * xxxx = new TVector(kcmaxhits*4);
2813 TVector & xxx = *xxxx;
2814 Int_t maxhits = kcmaxhits;
2817 AliTPChit * tpcHit=0;
2818 tpcHit = (AliTPChit*)FirstHit2(-1);
2819 Int_t currentIndex=0;
2820 Int_t lastrow=-1; //last writen row
2823 if (tpcHit==0) continue;
2824 sector=tpcHit->fSector; // sector number
2825 ipart=tpcHit->Track();
2829 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2830 Int_t index[3]={1,sector,0};
2831 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2832 if (currentrow<0) continue;
2833 if (lastrow<0) lastrow=currentrow;
2834 if (currentrow==lastrow){
2835 if ( currentIndex>=maxhits){
2837 xxx.ResizeTo(4*maxhits);
2839 xxx(currentIndex*4)=x[0];
2840 xxx(currentIndex*4+1)=x[1];
2841 xxx(currentIndex*4+2)=x[2];
2842 xxx(currentIndex*4+3)=tpcHit->fQ;
2846 if (currentIndex>2){
2858 for (Int_t index=0;index<currentIndex;index++){
2860 x=x2=x3=x4=xxx(index*4);
2868 sumy+=xxx(index*4+1);
2869 sumxy+=xxx(index*4+1)*x;
2870 sumx2y+=xxx(index*4+1)*x2;
2871 sumz+=xxx(index*4+2);
2872 sumxz+=xxx(index*4+2)*x;
2873 sumx2z+=xxx(index*4+2)*x2;
2874 sumq+=xxx(index*4+3);
2876 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2877 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2878 sumx2*(sumx*sumx3-sumx2*sumx2);
2880 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2881 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2882 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2883 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2885 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2886 sumx2*(sumx*sumx2y-sumx2*sumxy);
2887 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2888 sumx2*(sumx*sumx2z-sumx2*sumxz);
2890 Float_t y=detay/det+centralPad;
2891 Float_t z=detaz/det;
2892 Float_t by=detby/det; //y angle
2893 Float_t bz=detbz/det; //z angle
2894 sumy/=Float_t(currentIndex);
2895 sumz/=Float_t(currentIndex);
2897 AliComplexCluster cl;
2903 cl.fTracks[0]=ipart;
2905 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2906 if (row!=0) row->InsertCluster(&cl);
2909 } //end of calculating cluster for given row
2911 } // end of loop over hits
2915 //_______________________________________________________________________________
2916 void AliTPC::Digits2Reco(Int_t firstevent,Int_t lastevent)
2918 // produces rec points from digits and writes them on the root file
2919 // AliTPCclusters.root
2921 TDirectory *cwd = gDirectory;
2924 AliTPCParam *dig=(AliTPCParam *)gDirectory->Get("75x40_100x60");
2926 cout<<"AliTPC::Digits2Reco: TPC parameteres have been set"<<endl;
2928 if(!gSystem->Getenv("CONFIG_FILE")){
2929 out=TFile::Open("AliTPCclusters.root","recreate");
2935 tmp1=gSystem->Getenv("CONFIG_FILE_PREFIX");
2936 tmp2=gSystem->Getenv("CONFIG_OUTDIR");
2937 sprintf(tmp3,"%s%s/AliTPCclusters.root",tmp1,tmp2);
2938 out=TFile::Open(tmp3,"recreate");
2942 cout<<"AliTPC::Digits2Reco - determination of rec points begins"<<endl;
2945 for(Int_t iev=firstevent;iev<lastevent+1;iev++){
2947 printf("Processing event %d\n",iev);
2948 Digits2Clusters(out,iev);
2950 cout<<"AliTPC::Digits2Reco - determination of rec points ended"<<endl;