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.27 2001/01/13 17:29:33 kowal2
19 Sun compiler correction
21 Revision 1.26 2001/01/10 07:59:43 kowal2
22 Corrections to load points from the noncompressed hits.
24 Revision 1.25 2000/11/02 07:25:31 kowal2
25 Changes due to the new hit structure.
28 Revision 1.24 2000/10/05 16:06:09 kowal2
29 Forward declarations. Changes due to a new class AliComplexCluster.
31 Revision 1.23 2000/10/02 21:28:18 fca
32 Removal of useless dependecies via forward declarations
34 Revision 1.22 2000/07/10 20:57:39 hristov
35 Update of TPC code and macros by M.Kowalski
37 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
38 Changes to obey the coding rules
40 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
41 Splitted from AliTPCtracking
43 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
44 Changed parameter settings
46 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
48 Defaults loaded automatically (hard-wired)
49 Optional parameters can be set via macro called in the constructor
51 Revision 1.19 2000/04/18 19:00:59 fca
52 Small bug fixes to TPC files
54 Revision 1.18 2000/04/17 09:37:33 kowal2
55 removed obsolete AliTPCDigitsDisplay.C
57 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
59 New, experimental data structure from M. Ivanov
60 New tracking algorithm
61 Different pad geometry for different sectors
62 Digitization rewritten
64 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
65 Not used anymore - removed
67 Revision 1.17 2000/01/19 17:17:30 fca
68 Introducing a list of lists of hits -- more hits allowed for detector now
70 Revision 1.16 1999/11/05 09:29:23 fca
71 Accept only signals > 0
73 Revision 1.15 1999/10/08 06:26:53 fca
74 Removed ClustersIndex - not used anymore
76 Revision 1.14 1999/09/29 09:24:33 fca
77 Introduction of the Copyright and cvs Log
81 ///////////////////////////////////////////////////////////////////////////////
83 // Time Projection Chamber //
84 // This class contains the basic functions for the Time Projection Chamber //
85 // detector. Functions specific to one particular geometry are //
86 // contained in the derived classes //
90 <img src="picts/AliTPCClass.gif">
95 ///////////////////////////////////////////////////////////////////////////////
103 #include <TGeometry.h>
106 #include <TObjectTable.h>
107 #include "TParticle.h"
111 #include <iostream.h>
118 #include "AliTPCParamSR.h"
119 #include "AliTPCPRF2D.h"
120 #include "AliTPCRF1D.h"
121 #include "AliDigits.h"
122 #include "AliSimDigits.h"
123 #include "AliTPCTrackHits.h"
124 #include "AliPoints.h"
125 #include "AliArrayBranch.h"
128 #include "AliTPCDigitsArray.h"
129 #include "AliComplexCluster.h"
130 #include "AliClusters.h"
131 #include "AliTPCClustersRow.h"
132 #include "AliTPCClustersArray.h"
134 #include "AliTPCcluster.h"
135 #include "AliTPCclusterer.h"
136 #include "AliTPCtracker.h"
138 #include <TInterpreter.h>
145 //_____________________________________________________________________________
149 // Default constructor
164 //_____________________________________________________________________________
165 AliTPC::AliTPC(const char *name, const char *title)
166 : AliDetector(name,title)
169 // Standard constructor
173 // Initialise arrays of hits and digits
174 fHits = new TClonesArray("AliTPChit", 176);
175 gAlice->AddHitList(fHits);
180 fTrackHits = new AliTPCTrackHits; //MI - 13.09.2000
181 fTrackHits->SetHitPrecision(0.002);
182 fTrackHits->SetStepPrecision(0.003);
183 fTrackHits->SetMaxDistance(100);
187 // Initialise counters
193 // Initialise color attributes
194 SetMarkerColor(kYellow);
197 // Set TPC parameters
200 if (!strcmp(title,"Default")) {
201 fTPCParam = new AliTPCParamSR;
203 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
209 //_____________________________________________________________________________
220 delete fTrackHits; //MI 15.09.2000
223 //_____________________________________________________________________________
224 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
227 // Add a hit to the list
229 // TClonesArray &lhits = *fHits;
230 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
232 TClonesArray &lhits = *fHits;
233 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
236 AddHit2(track,vol,hits);
239 //_____________________________________________________________________________
240 void AliTPC::BuildGeometry()
244 // Build TPC ROOT TNode geometry for the event display
249 const int kColorTPC=19;
250 char name[5], title[25];
251 const Double_t kDegrad=TMath::Pi()/180;
252 const Double_t kRaddeg=180./TMath::Pi();
255 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
256 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
258 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
259 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
261 Int_t nLo = fTPCParam->GetNInnerSector()/2;
262 Int_t nHi = fTPCParam->GetNOuterSector()/2;
264 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
265 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
266 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
267 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
270 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
271 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
277 // Get ALICE top node
280 nTop=gAlice->GetGeometry()->GetNode("alice");
284 rl = fTPCParam->GetInnerRadiusLow();
285 ru = fTPCParam->GetInnerRadiusUp();
289 sprintf(name,"LS%2.2d",i);
291 sprintf(title,"TPC low sector %3d",i);
294 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
295 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
296 tubs->SetNumberOfDivisions(1);
298 nNode = new TNode(name,title,name,0,0,0,"");
299 nNode->SetLineColor(kColorTPC);
305 rl = fTPCParam->GetOuterRadiusLow();
306 ru = fTPCParam->GetOuterRadiusUp();
309 sprintf(name,"US%2.2d",i);
311 sprintf(title,"TPC upper sector %d",i);
313 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
314 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
315 tubs->SetNumberOfDivisions(1);
317 nNode = new TNode(name,title,name,0,0,0,"");
318 nNode->SetLineColor(kColorTPC);
324 //_____________________________________________________________________________
325 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
328 // Calculate distance from TPC to mouse on the display
334 void AliTPC::Clusters2Tracks(TFile *of) {
335 //-----------------------------------------------------------------
336 // This is a track finder.
337 //-----------------------------------------------------------------
338 AliTPCtracker::Clusters2Tracks(fTPCParam,of);
341 //_____________________________________________________________________________
342 void AliTPC::CreateMaterials()
344 //-----------------------------------------------
345 // Create Materials for for TPC simulations
346 //-----------------------------------------------
348 //-----------------------------------------------------------------
349 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
350 //-----------------------------------------------------------------
352 Int_t iSXFLD=gAlice->Field()->Integ();
353 Float_t sXMGMX=gAlice->Field()->Max();
355 Float_t amat[5]; // atomic numbers
356 Float_t zmat[5]; // z
357 Float_t wmat[5]; // proportions
363 //***************** Gases *************************
365 //-------------------------------------------------
367 //-------------------------------------------------
378 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
388 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
391 //--------------------------------------------------------------
393 //--------------------------------------------------------------
410 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
412 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
427 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
429 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
445 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
447 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
449 //----------------------------------------------------------------
450 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
451 //----------------------------------------------------------------
457 Float_t rho,absl,X0,buf[1];
461 for(nc = 0;nc<fNoComp;nc++)
464 // retrive material constants
466 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
471 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
473 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
474 density += fMixtProp[nc]*rho; // density of the mixture
478 // mixture proportions by weight!
480 for(nc = 0;nc<fNoComp;nc++)
483 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
485 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
486 apure[nnc] : amol[nnc])/am;
490 // Drift gases 1 - nonsensitive, 2 - sensitive
492 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
493 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
502 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
505 //----------------------------------------------------------------------
507 //----------------------------------------------------------------------
529 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
551 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
569 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
587 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
589 // G10 60% SiO2 + 40% epoxy, I use A and Z for SiO2
602 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz
604 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
606 AliMaterial(39,"G10",amat[0],zmat[0],density,999.,999.);
615 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
624 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
633 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
651 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
670 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
674 //----------------------------------------------------------
675 // tracking media for gases
676 //----------------------------------------------------------
678 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
679 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
680 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
681 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
683 //-----------------------------------------------------------
684 // tracking media for solids
685 //-----------------------------------------------------------
687 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
688 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
689 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
690 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
691 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
692 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
693 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
694 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
695 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
696 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
701 void AliTPC::Digits2Clusters(TFile *of)
703 //-----------------------------------------------------------------
704 // This is a simple cluster finder.
705 //-----------------------------------------------------------------
706 AliTPCclusterer::Digits2Clusters(fTPCParam,of);
709 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
710 extern Double_t SigmaZ2(Double_t, Double_t);
711 //_____________________________________________________________________________
712 void AliTPC::Hits2Clusters(TFile *of)
714 //--------------------------------------------------------
715 // TPC simple cluster generator from hits
716 // obtained from the TPC Fast Simulator
717 // The point errors are taken from the parametrization
718 //--------------------------------------------------------
720 //-----------------------------------------------------------------
721 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
722 //-----------------------------------------------------------------
723 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
724 // Jouri.Belikov@cern.ch
725 //----------------------------------------------------------------
727 /////////////////////////////////////////////////////////////////////////////
729 //---------------------------------------------------------------------
730 // ALICE TPC Cluster Parameters
731 //--------------------------------------------------------------------
735 // Cluster width in rphi
736 const Float_t kACrphi=0.18322;
737 const Float_t kBCrphi=0.59551e-3;
738 const Float_t kCCrphi=0.60952e-1;
739 // Cluster width in z
740 const Float_t kACz=0.19081;
741 const Float_t kBCz=0.55938e-3;
742 const Float_t kCCz=0.30428;
744 TDirectory *savedir=gDirectory;
747 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
752 printf("AliTPCParam MUST be created firstly\n");
756 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
758 TParticle *particle; // pointer to a given particle
759 AliTPChit *tpcHit; // pointer to a sigle TPC hit
763 Float_t pl,pt,tanth,rpad,ratio;
766 //---------------------------------------------------------------
767 // Get the access to the tracks
768 //---------------------------------------------------------------
770 TTree *tH = gAlice->TreeH();
771 Stat_t ntracks = tH->GetEntries();
773 //Switch to the output file
776 fTPCParam->Write(fTPCParam->GetTitle());
777 AliTPCClustersArray carray;
778 carray.Setup(fTPCParam);
779 carray.SetClusterType("AliTPCcluster");
782 Int_t nclusters=0; //cluster counter
784 //------------------------------------------------------------
785 // Loop over all sectors (72 sectors for 20 deg
786 // segmentation for both lower and upper sectors)
787 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
788 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
790 // First cluster for sector 0 starts at "0"
791 //------------------------------------------------------------
793 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
795 fTPCParam->AdjustCosSin(isec,cph,sph);
797 //------------------------------------------------------------
799 //------------------------------------------------------------
801 for(Int_t track=0;track<ntracks;track++){
805 // Get number of the TPC hits
807 // nhits=fHits->GetEntriesFast();
810 tpcHit = (AliTPChit*)FirstHit(-1);
814 // for(Int_t hit=0;hit<nhits;hit++){
815 //tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
819 if (tpcHit->fQ == 0.) {
820 tpcHit = (AliTPChit*) NextHit();
821 continue; //information about track (I.Belikov)
823 sector=tpcHit->fSector; // sector number
826 // if(sector != isec) continue; //terminate iteration
829 tpcHit = (AliTPChit*) NextHit();
832 ipart=tpcHit->Track();
833 particle=gAlice->Particle(ipart);
836 if(pt < 1.e-9) pt=1.e-9;
838 tanth = TMath::Abs(tanth);
839 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
840 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
842 // space-point resolutions
844 sigmaRphi=SigmaY2(rpad,tanth,pt);
845 sigmaZ =SigmaZ2(rpad,tanth );
849 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
850 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
852 // temporary protection
854 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
855 if(sigmaZ < 0.) sigmaZ=0.4e-3;
856 if(clRphi < 0.) clRphi=2.5e-3;
857 if(clZ < 0.) clZ=2.5e-5;
862 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
863 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
865 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
866 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
867 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
868 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
869 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
870 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
871 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
872 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
873 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
874 xyz[2]=sigmaRphi; // fSigmaY2
875 xyz[3]=sigmaZ; // fSigmaZ2
876 xyz[4]=tpcHit->fQ; // q
878 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
879 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
881 Int_t tracks[3]={tpcHit->Track(), -1, -1};
882 AliTPCcluster cluster(tracks,xyz);
884 clrow->InsertCluster(&cluster); nclusters++;
886 tpcHit = (AliTPChit*)NextHit();
889 } // end of loop over hits
891 } // end of loop over tracks
893 Int_t nrows=fTPCParam->GetNRow(isec);
894 for (Int_t irow=0; irow<nrows; irow++) {
895 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
896 if (!clrow) continue;
897 carray.StoreRow(isec,irow);
898 carray.ClearRow(isec,irow);
901 } // end of loop over sectors
903 cerr<<"Number of made clusters : "<<nclusters<<" \n";
905 carray.GetTree()->Write();
907 savedir->cd(); //switch back to the input file
911 //_________________________________________________________________
912 void AliTPC::Hits2ExactClustersSector(Int_t isec)
914 //--------------------------------------------------------
915 //calculate exact cross point of track and given pad row
916 //resulting values are expressed in "digit" coordinata
917 //--------------------------------------------------------
919 //-----------------------------------------------------------------
920 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
921 //-----------------------------------------------------------------
923 if (fClustersArray==0){
927 TParticle *particle; // pointer to a given particle
928 AliTPChit *tpcHit; // pointer to a sigle TPC hit
931 const Int_t kcmaxhits=30000;
932 TVector * xxxx = new TVector(kcmaxhits*4);
933 TVector & xxx = *xxxx;
934 Int_t maxhits = kcmaxhits;
935 //construct array for each padrow
936 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
937 fClustersArray->CreateRow(isec,i);
939 //---------------------------------------------------------------
940 // Get the access to the tracks
941 //---------------------------------------------------------------
943 TTree *tH = gAlice->TreeH();
944 Stat_t ntracks = tH->GetEntries();
945 Int_t npart = gAlice->GetNtrack();
947 //------------------------------------------------------------
949 //------------------------------------------------------------
951 for(Int_t track=0;track<ntracks;track++){
955 // Get number of the TPC hits and a pointer
958 nhits=fHits->GetEntriesFast();
962 Int_t currentIndex=0;
963 Int_t lastrow=-1; //last writen row
964 for(Int_t hit=0;hit<nhits;hit++){
965 tpcHit=(AliTPChit*)fHits->UncheckedAt(hit);
966 if (tpcHit==0) continue;
967 sector=tpcHit->fSector; // sector number
968 if(sector != isec) continue;
969 ipart=tpcHit->Track();
970 if (ipart<npart) particle=gAlice->Particle(ipart);
974 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
975 Int_t index[3]={1,isec,0};
976 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
977 if (currentrow<0) continue;
978 if (lastrow<0) lastrow=currentrow;
979 if (currentrow==lastrow){
980 if ( currentIndex>=maxhits){
982 xxx.ResizeTo(4*maxhits);
984 xxx(currentIndex*4)=x[0];
985 xxx(currentIndex*4+1)=x[1];
986 xxx(currentIndex*4+2)=x[2];
987 xxx(currentIndex*4+3)=tpcHit->fQ;
1003 for (Int_t index=0;index<currentIndex;index++){
1005 x=x2=x3=x4=xxx(index*4);
1013 sumy+=xxx(index*4+1);
1014 sumxy+=xxx(index*4+1)*x;
1015 sumx2y+=xxx(index*4+1)*x2;
1016 sumz+=xxx(index*4+2);
1017 sumxz+=xxx(index*4+2)*x;
1018 sumx2z+=xxx(index*4+2)*x2;
1019 sumq+=xxx(index*4+3);
1021 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1022 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1023 sumx2*(sumx*sumx3-sumx2*sumx2);
1025 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1026 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1027 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1028 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1030 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1031 sumx2*(sumx*sumx2y-sumx2*sumxy);
1032 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1033 sumx2*(sumx*sumx2z-sumx2*sumxz);
1035 Float_t y=detay/det+centralPad;
1036 Float_t z=detaz/det;
1037 Float_t by=detby/det; //y angle
1038 Float_t bz=detbz/det; //z angle
1039 sumy/=Float_t(currentIndex);
1040 sumz/=Float_t(currentIndex);
1041 AliComplexCluster cl;
1047 cl.fTracks[0]=ipart;
1049 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1050 if (row!=0) row->InsertCluster(&cl);
1053 } //end of calculating cluster for given row
1057 } // end of loop over hits
1058 } // end of loop over tracks
1059 //write padrows to tree
1060 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1061 fClustersArray->StoreRow(isec,ii);
1062 fClustersArray->ClearRow(isec,ii);
1067 //___________________________________________
1068 void AliTPC::SDigits2Digits()
1070 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1071 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1072 AliTPCPRF2D * prfouter = new AliTPCPRF2D;
1073 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1075 TDirectory *cwd = gDirectory;
1076 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1077 rf->SetOffset(3*param->GetZSigma());
1079 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1081 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n";
1084 prfinner->Read("prf_07504_Gati_056068_d02");
1085 prfouter->Read("prf_10006_Gati_047051_d03");
1089 param->SetInnerPRF(prfinner);
1090 param->SetOuterPRF(prfouter);
1091 param->SetTimeRF(rf);
1095 cerr<<"Digitizing TPC...\n";
1097 //setup TPCDigitsArray
1098 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1099 arr->SetClass("AliSimDigits");
1103 arr->MakeTree(fDigitsFile);
1105 SetDigitsArray(arr);
1109 // Hits2DigitsSector(1);
1110 // Hits2DigitsSector(2);
1111 // Hits2DigitsSector(3);
1112 // Hits2DigitsSector(1+18);
1113 // Hits2DigitsSector(2+18);
1114 // Hits2DigitsSector(3+18);
1116 // Hits2DigitsSector(36+1);
1117 // Hits2DigitsSector(36+2);
1118 // Hits2DigitsSector(36+3);
1119 // Hits2DigitsSector(36+1+18);
1120 // Hits2DigitsSector(36+2+18);
1121 // Hits2DigitsSector(36+3+18);
1126 sprintf(treeName,"TreeD_%s",param->GetTitle());
1127 GetDigitsArray()->GetTree()->Write(treeName,TObject::kOverwrite);
1130 //__________________________________________________________________
1131 void AliTPC::Hits2Digits()
1133 //----------------------------------------------------
1134 // Loop over all sectors
1135 //----------------------------------------------------
1138 printf("AliTPCParam MUST be created firstly\n");
1142 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) Hits2DigitsSector(isec);
1147 //_____________________________________________________________________________
1148 void AliTPC::Hits2DigitsSector(Int_t isec)
1150 //-------------------------------------------------------------------
1151 // TPC conversion from hits to digits.
1152 //-------------------------------------------------------------------
1154 //-----------------------------------------------------------------
1155 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1156 //-----------------------------------------------------------------
1158 //-------------------------------------------------------
1159 // Get the access to the track hits
1160 //-------------------------------------------------------
1163 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1164 Stat_t ntracks = tH->GetEntries();
1168 //-------------------------------------------
1169 // Only if there are any tracks...
1170 //-------------------------------------------
1174 printf("*** Processing sector number %d ***\n",isec);
1176 Int_t nrows =fTPCParam->GetNRow(isec);
1178 row= new TObjArray* [nrows];
1180 MakeSector(isec,nrows,tH,ntracks,row);
1182 //--------------------------------------------------------
1183 // Digitize this sector, row by row
1184 // row[i] is the pointer to the TObjArray of TVectors,
1185 // each one containing electrons accepted on this
1186 // row, assigned into tracks
1187 //--------------------------------------------------------
1191 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1193 for (i=0;i<nrows;i++){
1195 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1197 DigitizeRow(i,isec,row);
1199 fDigitsArray->StoreRow(isec,i);
1201 Int_t ndig = dig->GetDigitSize();
1203 printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1205 fDigitsArray->ClearRow(isec,i);
1208 } // end of the sector digitization
1210 for(i=0;i<nrows;i++){
1215 delete [] row; // delete the array of pointers to TObjArray-s
1219 } // end of Hits2DigitsSector
1222 //_____________________________________________________________________________
1223 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1225 //-----------------------------------------------------------
1226 // Single row digitization, coupling from the neighbouring
1227 // rows taken into account
1228 //-----------------------------------------------------------
1230 //-----------------------------------------------------------------
1231 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1232 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1233 //-----------------------------------------------------------------
1236 Float_t zerosup = fTPCParam->GetZeroSup();
1237 Int_t nrows =fTPCParam->GetNRow(isec);
1238 fCurrentIndex[1]= isec;
1241 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1242 Int_t nofTbins = fTPCParam->GetMaxTBin();
1243 Int_t indexRange[4];
1245 // Integrated signal for this row
1246 // and a single track signal
1248 TMatrix *m1 = new TMatrix(0,nofPads,0,nofTbins); // integrated
1249 TMatrix *m2 = new TMatrix(0,nofPads,0,nofTbins); // single
1251 TMatrix &total = *m1;
1253 // Array of pointers to the label-signal list
1255 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1256 Float_t **pList = new Float_t* [nofDigits];
1260 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1264 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1265 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1266 for (Int_t row= row1;row<=row2;row++){
1267 Int_t nTracks= rows[row]->GetEntries();
1268 for (i1=0;i1<nTracks;i1++){
1269 fCurrentIndex[2]= row;
1270 fCurrentIndex[3]=irow;
1272 m2->Zero(); // clear single track signal matrix
1273 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1274 GetList(trackLabel,nofPads,m2,indexRange,pList);
1276 else GetSignal(rows[row],i1,0,m1,indexRange);
1282 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1283 for(Int_t ip=0;ip<nofPads;ip++){
1284 for(Int_t it=0;it<nofTbins;it++){
1286 Float_t q = total(ip,it);
1288 Int_t gi =it*nofPads+ip; // global index
1290 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
1294 if(q <=zerosup) continue; // do not fill zeros
1295 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1298 // "real" signal or electronic noise (list = -1)?
1301 for(Int_t j1=0;j1<3;j1++){
1302 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -1;
1307 <A NAME="AliDigits"></A>
1308 using of AliDigits object
1311 dig->SetDigitFast((Short_t)q,it,ip);
1312 if (fDigitsArray->IsSimulated())
1314 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1315 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1316 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1320 } // end of loop over time buckets
1321 } // end of lop over pads
1324 // This row has been digitized, delete nonused stuff
1327 for(lp=0;lp<nofDigits;lp++){
1328 if(pList[lp]) delete [] pList[lp];
1337 } // end of DigitizeRow
1339 //_____________________________________________________________________________
1341 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr, TMatrix *m1, TMatrix *m2,
1345 //---------------------------------------------------------------
1346 // Calculates 2-D signal (pad,time) for a single track,
1347 // returns a pointer to the signal matrix and the track label
1348 // No digitization is performed at this level!!!
1349 //---------------------------------------------------------------
1351 //-----------------------------------------------------------------
1352 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1353 // Modified: Marian Ivanov
1354 //-----------------------------------------------------------------
1358 tv = (TVector*)p1->At(ntr); // pointer to a track
1361 Float_t label = v(0);
1362 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1364 Int_t nElectrons = (tv->GetNrows()-1)/4;
1365 indexRange[0]=9999; // min pad
1366 indexRange[1]=-1; // max pad
1367 indexRange[2]=9999; //min time
1368 indexRange[3]=-1; // max time
1370 // Float_t IneffFactor = 0.5; // inefficiency in the gain close to the edge, as above
1372 TMatrix &signal = *m1;
1373 TMatrix &total = *m2;
1375 // Loop over all electrons
1377 for(Int_t nel=0; nel<nElectrons; nel++){
1379 Float_t aval = v(idx+4);
1380 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1381 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1382 Int_t n = fTPCParam->CalcResponse(xyz,fCurrentIndex,fCurrentIndex[3]);
1384 if (n>0) for (Int_t i =0; i<n; i++){
1385 Int_t *index = fTPCParam->GetResBin(i);
1386 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1387 if ( ( pad<(fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]))) && (pad>0)) {
1388 Int_t time=index[2];
1389 Float_t weight = fTPCParam->GetResWeight(i); //we normalise response to ADC channel
1390 weight *= eltoadcfac;
1392 if (m1!=0) signal(pad,time)+=weight;
1393 total(pad,time)+=weight;
1394 indexRange[0]=TMath::Min(indexRange[0],pad);
1395 indexRange[1]=TMath::Max(indexRange[1],pad);
1396 indexRange[2]=TMath::Min(indexRange[2],time);
1397 indexRange[3]=TMath::Max(indexRange[3],time);
1400 } // end of loop over electrons
1402 return label; // returns track label when finished
1405 //_____________________________________________________________________________
1406 void AliTPC::GetList(Float_t label,Int_t np,TMatrix *m,Int_t *indexRange,
1409 //----------------------------------------------------------------------
1410 // Updates the list of tracks contributing to digits for a given row
1411 //----------------------------------------------------------------------
1413 //-----------------------------------------------------------------
1414 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1415 //-----------------------------------------------------------------
1417 TMatrix &signal = *m;
1419 // lop over nonzero digits
1421 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1422 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1425 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1427 if(signal(ip,it)<0.5) continue;
1430 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1432 if(!pList[globalIndex]){
1435 // Create new list (6 elements - 3 signals and 3 labels),
1438 pList[globalIndex] = new Float_t [6];
1442 *pList[globalIndex] = -1.;
1443 *(pList[globalIndex]+1) = -1.;
1444 *(pList[globalIndex]+2) = -1.;
1445 *(pList[globalIndex]+3) = -1.;
1446 *(pList[globalIndex]+4) = -1.;
1447 *(pList[globalIndex]+5) = -1.;
1450 *pList[globalIndex] = label;
1451 *(pList[globalIndex]+3) = signal(ip,it);
1455 // check the signal magnitude
1457 Float_t highest = *(pList[globalIndex]+3);
1458 Float_t middle = *(pList[globalIndex]+4);
1459 Float_t lowest = *(pList[globalIndex]+5);
1462 // compare the new signal with already existing list
1465 if(signal(ip,it)<lowest) continue; // neglect this track
1469 if (signal(ip,it)>highest){
1470 *(pList[globalIndex]+5) = middle;
1471 *(pList[globalIndex]+4) = highest;
1472 *(pList[globalIndex]+3) = signal(ip,it);
1474 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1475 *(pList[globalIndex]+1) = *pList[globalIndex];
1476 *pList[globalIndex] = label;
1478 else if (signal(ip,it)>middle){
1479 *(pList[globalIndex]+5) = middle;
1480 *(pList[globalIndex]+4) = signal(ip,it);
1482 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1483 *(pList[globalIndex]+1) = label;
1486 *(pList[globalIndex]+5) = signal(ip,it);
1487 *(pList[globalIndex]+2) = label;
1491 } // end of loop over pads
1492 } // end of loop over time bins
1497 //___________________________________________________________________
1498 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1499 Stat_t ntracks,TObjArray **row)
1502 //-----------------------------------------------------------------
1503 // Prepares the sector digitization, creates the vectors of
1504 // tracks for each row of this sector. The track vector
1505 // contains the track label and the position of electrons.
1506 //-----------------------------------------------------------------
1508 //-----------------------------------------------------------------
1509 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1510 //-----------------------------------------------------------------
1512 Float_t gasgain = fTPCParam->GetGasGain();
1516 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1519 if (fHitType&2) branch = TH->GetBranch("TPC2");
1520 else branch = TH->GetBranch("TPC");
1523 //----------------------------------------------
1524 // Create TObjArray-s, one for each row,
1525 // each TObjArray will store the TVectors
1526 // of electrons, one TVector per each track.
1527 //----------------------------------------------
1529 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
1530 TVector **tracks = new TVector* [nrows]; //pointers to the track vectors
1531 for(i=0; i<nrows; i++){
1532 row[i] = new TObjArray;
1539 //--------------------------------------------------------------------
1540 // Loop over tracks, the "track" contains the full history
1541 //--------------------------------------------------------------------
1543 Int_t previousTrack,currentTrack;
1544 previousTrack = -1; // nothing to store so far!
1546 for(Int_t track=0;track<ntracks;track++){
1547 Bool_t isInSector=kTRUE;
1552 TBranch * br = TH->GetBranch("fTrackHitsInfo");
1553 br->GetEvent(track);
1554 AliObjectArray * ar = fTrackHits->fTrackHitsInfo;
1555 for (UInt_t j=0;j<ar->GetSize();j++){
1556 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==isec) isInSector=kTRUE;
1559 if (!isInSector) continue;
1561 branch->GetEntry(track); // get next track
1565 tpcHit = (AliTPChit*)FirstHit(-1);
1567 //--------------------------------------------------------------
1569 //--------------------------------------------------------------
1574 Int_t sector=tpcHit->fSector; // sector number
1575 // if(sector != isec) continue;
1577 tpcHit = (AliTPChit*) NextHit();
1581 currentTrack = tpcHit->Track(); // track number
1584 if(currentTrack != previousTrack){
1586 // store already filled fTrack
1588 for(i=0;i<nrows;i++){
1589 if(previousTrack != -1){
1590 if(nofElectrons[i]>0){
1591 TVector &v = *tracks[i];
1592 v(0) = previousTrack;
1593 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1594 row[i]->Add(tracks[i]);
1597 delete tracks[i]; // delete empty TVector
1603 tracks[i] = new TVector(481); // TVectors for the next fTrack
1605 } // end of loop over rows
1607 previousTrack=currentTrack; // update track label
1610 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
1612 //---------------------------------------------------
1613 // Calculate the electron attachment probability
1614 //---------------------------------------------------
1617 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
1618 /fTPCParam->GetDriftV();
1620 Float_t attProb = fTPCParam->GetAttCoef()*
1621 fTPCParam->GetOxyCont()*time; // fraction!
1623 //-----------------------------------------------
1624 // Loop over electrons
1625 //-----------------------------------------------
1628 for(Int_t nel=0;nel<qI;nel++){
1629 // skip if electron lost due to the attachment
1630 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
1634 xyz[3]= (Float_t) (-gasgain*TMath::Log(gRandom->Rndm()));
1637 TransportElectron(xyz,index); //MI change -august
1639 fTPCParam->GetPadRow(xyz,index); //MI change august
1640 rowNumber = index[2];
1641 //transform position to local digit coordinates
1642 //relative to nearest pad row
1643 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
1644 nofElectrons[rowNumber]++;
1645 //----------------------------------
1646 // Expand vector if necessary
1647 //----------------------------------
1648 if(nofElectrons[rowNumber]>120){
1649 Int_t range = tracks[rowNumber]->GetNrows();
1650 if((nofElectrons[rowNumber])>(range-1)/4){
1652 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
1656 TVector &v = *tracks[rowNumber];
1657 Int_t idx = 4*nofElectrons[rowNumber]-3;
1659 v(idx)= xyz[0]; // X - pad row coordinate
1660 v(idx+1)=xyz[1]; // Y - pad coordinate (along the pad-row)
1661 v(idx+2)=xyz[2]; // Z - time bin coordinate
1662 v(idx+3)=xyz[3]; // avalanche size
1663 } // end of loop over electrons
1665 tpcHit = (AliTPChit*)NextHit();
1667 } // end of loop over hits
1668 } // end of loop over tracks
1671 // store remaining track (the last one) if not empty
1674 for(i=0;i<nrows;i++){
1675 if(nofElectrons[i]>0){
1676 TVector &v = *tracks[i];
1677 v(0) = previousTrack;
1678 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
1679 row[i]->Add(tracks[i]);
1688 delete [] nofElectrons;
1691 } // end of MakeSector
1694 //_____________________________________________________________________________
1698 // Initialise TPC detector after definition of geometry
1703 for(i=0;i<35;i++) printf("*");
1704 printf(" TPC_INIT ");
1705 for(i=0;i<35;i++) printf("*");
1708 for(i=0;i<80;i++) printf("*");
1712 //_____________________________________________________________________________
1713 void AliTPC::MakeBranch(Option_t* option, char *file)
1716 // Create Tree branches for the TPC.
1718 Int_t buffersize = 4000;
1719 char branchname[10];
1720 sprintf(branchname,"%s",GetName());
1722 AliDetector::MakeBranch(option,file);
1724 char *d = strstr(option,"D");
1726 if (fDigits && gAlice->TreeD() && d) {
1727 gAlice->MakeBranchInTree(gAlice->TreeD(),
1728 branchname, &fDigits, buffersize, file) ;
1731 if (fHitType&2) MakeBranch2(option,file); // MI change 14.09.2000
1734 //_____________________________________________________________________________
1735 void AliTPC::ResetDigits()
1738 // Reset number of digits and the digits array for this detector
1741 if (fDigits) fDigits->Clear();
1744 //_____________________________________________________________________________
1745 void AliTPC::SetSecAL(Int_t sec)
1747 //---------------------------------------------------
1748 // Activate/deactivate selection for lower sectors
1749 //---------------------------------------------------
1751 //-----------------------------------------------------------------
1752 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1753 //-----------------------------------------------------------------
1758 //_____________________________________________________________________________
1759 void AliTPC::SetSecAU(Int_t sec)
1761 //----------------------------------------------------
1762 // Activate/deactivate selection for upper sectors
1763 //---------------------------------------------------
1765 //-----------------------------------------------------------------
1766 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1767 //-----------------------------------------------------------------
1772 //_____________________________________________________________________________
1773 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
1775 //----------------------------------------
1776 // Select active lower sectors
1777 //----------------------------------------
1779 //-----------------------------------------------------------------
1780 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1781 //-----------------------------------------------------------------
1791 //_____________________________________________________________________________
1792 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
1793 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
1794 Int_t s11 , Int_t s12)
1796 //--------------------------------
1797 // Select active upper sectors
1798 //--------------------------------
1800 //-----------------------------------------------------------------
1801 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1802 //-----------------------------------------------------------------
1818 //_____________________________________________________________________________
1819 void AliTPC::SetSens(Int_t sens)
1822 //-------------------------------------------------------------
1823 // Activates/deactivates the sensitive strips at the center of
1824 // the pad row -- this is for the space-point resolution calculations
1825 //-------------------------------------------------------------
1827 //-----------------------------------------------------------------
1828 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1829 //-----------------------------------------------------------------
1835 void AliTPC::SetSide(Float_t side=0.)
1837 // choice of the TPC side
1842 //____________________________________________________________________________
1843 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
1844 Float_t p2,Float_t p3)
1847 // gax mixture definition
1861 //_____________________________________________________________________________
1863 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
1866 // electron transport taking into account:
1868 // 2.ExB at the wires
1869 // 3. nonisochronity
1871 // xyz and index must be already transformed to system 1
1874 fTPCParam->Transform1to2(xyz,index);
1877 Float_t driftl=xyz[2];
1878 if(driftl<0.01) driftl=0.01;
1879 driftl=TMath::Sqrt(driftl);
1880 Float_t sigT = driftl*(fTPCParam->GetDiffT());
1881 Float_t sigL = driftl*(fTPCParam->GetDiffL());
1882 xyz[0]=gRandom->Gaus(xyz[0],sigT);
1883 xyz[1]=gRandom->Gaus(xyz[1],sigT);
1884 xyz[2]=gRandom->Gaus(xyz[2],sigL);
1888 if (fTPCParam->GetMWPCReadout()==kTRUE){
1890 fTPCParam->Transform2to2NearestWire(xyz,index);
1891 Float_t dx=xyz[0]-x1;
1892 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
1894 //add nonisochronity (not implemented yet)
1898 ClassImp(AliTPCdigit)
1900 //_____________________________________________________________________________
1901 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
1905 // Creates a TPC digit object
1907 fSector = digits[0];
1908 fPadRow = digits[1];
1911 fSignal = digits[4];
1917 //_____________________________________________________________________________
1918 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
1922 // Creates a TPC hit object
1933 //________________________________________________________________________
1934 // Additional code because of the AliTPCTrackHits
1936 void AliTPC::MakeBranch2(Option_t *option,char *file)
1939 // Create a new branch in the current Root Tree
1940 // The branch of fHits is automatically split
1941 // MI change 14.09.2000
1942 if (fHitType&2==0) return;
1943 char branchname[10];
1944 sprintf(branchname,"%s2",GetName());
1946 // Get the pointer to the header
1947 char *cH = strstr(option,"H");
1949 if (fTrackHits && gAlice->TreeH() && cH) {
1950 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHits,
1951 gAlice->TreeH(),fBufferSize,1);
1952 gAlice->TreeH()->GetListOfBranches()->Add(branch);
1953 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
1955 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
1956 TDirectory *wd = gDirectory;
1958 TIter next( b->GetListOfBranches());
1959 while ((b=(TBranch*)next())) {
1963 cout << "Diverting branch " << branchname << " to file " << file << endl;
1968 void AliTPC::SetTreeAddress()
1970 if (fHitType&1) AliDetector::SetTreeAddress();
1971 if (fHitType&2) SetTreeAddress2();
1974 void AliTPC::SetTreeAddress2()
1977 // Set branch address for the TrackHits Tree
1980 char branchname[20];
1981 sprintf(branchname,"%s2",GetName());
1983 // Branch address for hit tree
1984 TTree *treeH = gAlice->TreeH();
1986 branch = treeH->GetBranch(branchname);
1987 if (branch) branch->SetAddress(&fTrackHits);
1991 void AliTPC::FinishPrimary()
1993 if (fTrackHits) fTrackHits->FlushHitStack();
1997 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2000 // add hit to the list
2003 int primary = gAlice->GetPrimary(track);
2004 gAlice->Particle(primary)->SetBit(kKeepBit);
2008 gAlice->FlagTrack(track);
2010 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2011 //if (hit->fTrack!=rtrack)
2012 // cout<<"bad track number\n";
2014 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2015 hits[1],hits[2],(Int_t)hits[3]);
2018 void AliTPC::ResetHits()
2020 if (fHitType&1) AliDetector::ResetHits();
2021 if (fHitType&2) ResetHits2();
2024 void AliTPC::ResetHits2()
2028 if (fTrackHits) fTrackHits->Clear();
2031 AliHit* AliTPC::FirstHit(Int_t track)
2033 if (fHitType&2) return FirstHit2(track);
2034 return AliDetector::FirstHit(track);
2036 AliHit* AliTPC::NextHit()
2038 if (fHitType&2) return NextHit2();
2039 return AliDetector::NextHit();
2042 AliHit* AliTPC::FirstHit2(Int_t track)
2045 // Initialise the hit iterator
2046 // Return the address of the first hit for track
2047 // If track>=0 the track is read from disk
2048 // while if track<0 the first hit of the current
2049 // track is returned
2052 gAlice->ResetHits();
2053 gAlice->TreeH()->GetEvent(track);
2057 fTrackHits->First();
2058 return fTrackHits->GetHit();
2063 AliHit* AliTPC::NextHit2()
2066 //Return the next hit for the current track
2070 return fTrackHits->GetHit();
2076 void AliTPC::LoadPoints(Int_t)
2080 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2083 if(fHitType==1) AliDetector::LoadPoints(a);
2084 else LoadPoints2(a);
2091 void AliTPC::RemapTrackHitIDs(Int_t *map)
2093 if (!fTrackHits) return;
2094 AliObjectArray * arr = fTrackHits->fTrackHitsInfo;
2095 for (UInt_t i=0;i<arr->GetSize();i++){
2096 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2097 info->fTrackID = map[info->fTrackID];
2103 //_____________________________________________________________________________
2104 void AliTPC::LoadPoints2(Int_t)
2107 // Store x, y, z of all hits in memory
2109 if (fTrackHits == 0) return;
2111 Int_t nhits = fTrackHits->GetEntriesFast();
2112 if (nhits == 0) return;
2113 Int_t tracks = gAlice->GetNtrack();
2114 if (fPoints == 0) fPoints = new TObjArray(tracks);
2117 Int_t *ntrk=new Int_t[tracks];
2118 Int_t *limi=new Int_t[tracks];
2119 Float_t **coor=new Float_t*[tracks];
2120 for(Int_t i=0;i<tracks;i++) {
2126 AliPoints *points = 0;
2129 Int_t chunk=nhits/4+1;
2131 // Loop over all the hits and store their position
2133 ahit = FirstHit2(-1);
2134 //for (Int_t hit=0;hit<nhits;hit++) {
2136 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2137 trk=ahit->GetTrack();
2138 if(ntrk[trk]==limi[trk]) {
2140 // Initialise a new track
2141 fp=new Float_t[3*(limi[trk]+chunk)];
2143 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2144 delete [] coor[trk];
2151 fp[3*ntrk[trk] ] = ahit->X();
2152 fp[3*ntrk[trk]+1] = ahit->Y();
2153 fp[3*ntrk[trk]+2] = ahit->Z();
2158 for(trk=0; trk<tracks; ++trk) {
2160 points = new AliPoints();
2161 points->SetMarkerColor(GetMarkerColor());
2162 points->SetMarkerSize(GetMarkerSize());
2163 points->SetDetector(this);
2164 points->SetParticle(trk);
2165 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2166 fPoints->AddAt(points,trk);
2167 delete [] coor[trk];
2177 //_____________________________________________________________________________
2178 void AliTPC::LoadPoints3(Int_t)
2181 // Store x, y, z of all hits in memory
2182 // - only intersection point with pad row
2183 if (fTrackHits == 0) return;
2185 Int_t nhits = fTrackHits->GetEntriesFast();
2186 if (nhits == 0) return;
2187 Int_t tracks = gAlice->GetNtrack();
2188 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2189 fPoints->Expand(2*tracks);
2192 Int_t *ntrk=new Int_t[tracks];
2193 Int_t *limi=new Int_t[tracks];
2194 Float_t **coor=new Float_t*[tracks];
2195 for(Int_t i=0;i<tracks;i++) {
2201 AliPoints *points = 0;
2204 Int_t chunk=nhits/4+1;
2206 // Loop over all the hits and store their position
2208 ahit = FirstHit2(-1);
2209 //for (Int_t hit=0;hit<nhits;hit++) {
2213 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2214 trk=ahit->GetTrack();
2215 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2216 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2217 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2218 if (currentrow!=lastrow){
2219 lastrow = currentrow;
2220 //later calculate intersection point
2221 if(ntrk[trk]==limi[trk]) {
2223 // Initialise a new track
2224 fp=new Float_t[3*(limi[trk]+chunk)];
2226 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2227 delete [] coor[trk];
2234 fp[3*ntrk[trk] ] = ahit->X();
2235 fp[3*ntrk[trk]+1] = ahit->Y();
2236 fp[3*ntrk[trk]+2] = ahit->Z();
2243 for(trk=0; trk<tracks; ++trk) {
2245 points = new AliPoints();
2246 points->SetMarkerColor(GetMarkerColor()+1);
2247 points->SetMarkerStyle(5);
2248 points->SetMarkerSize(0.2);
2249 points->SetDetector(this);
2250 points->SetParticle(trk);
2251 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2252 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2253 fPoints->AddAt(points,tracks+trk);
2254 delete [] coor[trk];
2265 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
2269 //fill clones array with intersection of current point with the
2274 const Int_t kcmaxhits=30000;
2275 TVector * xxxx = new TVector(kcmaxhits*4);
2276 TVector & xxx = *xxxx;
2277 Int_t maxhits = kcmaxhits;
2280 AliTPChit * tpcHit=0;
2281 tpcHit = (AliTPChit*)FirstHit2(-1);
2282 Int_t currentIndex=0;
2283 Int_t lastrow=-1; //last writen row
2286 if (tpcHit==0) continue;
2287 sector=tpcHit->fSector; // sector number
2288 ipart=tpcHit->Track();
2292 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2293 Int_t index[3]={1,sector,0};
2294 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2295 if (currentrow<0) continue;
2296 if (lastrow<0) lastrow=currentrow;
2297 if (currentrow==lastrow){
2298 if ( currentIndex>=maxhits){
2300 xxx.ResizeTo(4*maxhits);
2302 xxx(currentIndex*4)=x[0];
2303 xxx(currentIndex*4+1)=x[1];
2304 xxx(currentIndex*4+2)=x[2];
2305 xxx(currentIndex*4+3)=tpcHit->fQ;
2309 if (currentIndex>2){
2321 for (Int_t index=0;index<currentIndex;index++){
2323 x=x2=x3=x4=xxx(index*4);
2331 sumy+=xxx(index*4+1);
2332 sumxy+=xxx(index*4+1)*x;
2333 sumx2y+=xxx(index*4+1)*x2;
2334 sumz+=xxx(index*4+2);
2335 sumxz+=xxx(index*4+2)*x;
2336 sumx2z+=xxx(index*4+2)*x2;
2337 sumq+=xxx(index*4+3);
2339 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2340 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2341 sumx2*(sumx*sumx3-sumx2*sumx2);
2343 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2344 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2345 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2346 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2348 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2349 sumx2*(sumx*sumx2y-sumx2*sumxy);
2350 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2351 sumx2*(sumx*sumx2z-sumx2*sumxz);
2353 Float_t y=detay/det+centralPad;
2354 Float_t z=detaz/det;
2355 Float_t by=detby/det; //y angle
2356 Float_t bz=detbz/det; //z angle
2357 sumy/=Float_t(currentIndex);
2358 sumz/=Float_t(currentIndex);
2360 AliComplexCluster cl;
2366 cl.fTracks[0]=ipart;
2368 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2369 if (row!=0) row->InsertCluster(&cl);
2372 } //end of calculating cluster for given row
2374 } // end of loop over hits