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 ///////////////////////////////////////////////////////////////////////////////
20 // Time Projection Chamber //
21 // This class contains the basic functions for the Time Projection Chamber //
22 // detector. Functions specific to one particular geometry are //
23 // contained in the derived classes //
27 <img src="picts/AliTPCClass.gif">
32 ///////////////////////////////////////////////////////////////////////////////
36 #include <Riostream.h>
40 #include <TGeometry.h>
41 #include <TInterpreter.h>
45 #include <TObjectTable.h>
46 #include <TParticle.h>
53 #include <TVirtualMC.h>
57 #include "AliArrayBranch.h"
58 #include "AliClusters.h"
59 #include "AliComplexCluster.h"
60 #include "AliDigits.h"
62 #include "AliPoints.h"
64 #include "AliRunLoader.h"
65 #include "AliSimDigits.h"
68 #include "AliTPCClustersArray.h"
69 #include "AliTPCClustersRow.h"
70 #include "AliTPCDigitsArray.h"
71 #include "AliTPCLoader.h"
72 #include "AliTPCPRF2D.h"
73 #include "AliTPCParamSR.h"
74 #include "AliTPCRF1D.h"
75 #include "AliTPCTrackHits.h"
76 #include "AliTPCTrackHitsV2.h"
77 #include "AliTPCcluster.h"
78 #include "AliTrackReference.h"
80 #include "AliTPCDigitizer.h"
85 //_____________________________________________________________________________
86 // helper class for fast matrix and vector manipulation - no range checking
87 // origin - Marian Ivanov
89 class AliTPCFastMatrix : public TMatrix {
91 AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb);
92 inline Float_t & UncheckedAt(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
93 inline Float_t UncheckedAtFast(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
96 AliTPCFastMatrix::AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb):
97 TMatrix(row_lwb, row_upb,col_lwb,col_upb)
100 //_____________________________________________________________________________
101 class AliTPCFastVector : public TVector {
103 AliTPCFastVector(Int_t size);
104 virtual ~AliTPCFastVector(){;}
105 inline Float_t & UncheckedAt(Int_t index) const {return fElements[index];} //fast acces
109 AliTPCFastVector::AliTPCFastVector(Int_t size):TVector(size){
112 //_____________________________________________________________________________
116 // Default constructor
127 fHitType = 2; //default CONTAINERS - based on ROOT structure
134 //_____________________________________________________________________________
135 AliTPC::AliTPC(const char *name, const char *title)
136 : AliDetector(name,title)
139 // Standard constructor
143 // Initialise arrays of hits and digits
144 fHits = new TClonesArray("AliTPChit", 176);
145 gAlice->GetMCApp()->AddHitList(fHits);
150 fTrackHits = new AliTPCTrackHitsV2;
151 fTrackHits->SetHitPrecision(0.002);
152 fTrackHits->SetStepPrecision(0.003);
153 fTrackHits->SetMaxDistance(100);
155 fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
156 fTrackHitsOld->SetHitPrecision(0.002);
157 fTrackHitsOld->SetStepPrecision(0.003);
158 fTrackHitsOld->SetMaxDistance(100);
165 // Initialise counters
171 // Initialise color attributes
172 SetMarkerColor(kYellow);
175 // Set TPC parameters
179 if (!strcmp(title,"Default")) {
180 fTPCParam = new AliTPCParamSR;
182 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
188 //_____________________________________________________________________________
199 delete fTrackHits; //MI 15.09.2000
200 delete fTrackHitsOld; //MI 10.12.2001
201 if (fNoiseTable) delete [] fNoiseTable;
205 //_____________________________________________________________________________
206 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
209 // Add a hit to the list
211 // TClonesArray &lhits = *fHits;
212 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
214 TClonesArray &lhits = *fHits;
215 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
218 AddHit2(track,vol,hits);
221 //_____________________________________________________________________________
222 void AliTPC::BuildGeometry()
226 // Build TPC ROOT TNode geometry for the event display
231 const int kColorTPC=19;
232 char name[5], title[25];
233 const Double_t kDegrad=TMath::Pi()/180;
234 const Double_t kRaddeg=180./TMath::Pi();
237 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
238 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
240 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
241 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
243 Int_t nLo = fTPCParam->GetNInnerSector()/2;
244 Int_t nHi = fTPCParam->GetNOuterSector()/2;
246 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
247 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
248 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
249 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
252 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
253 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
259 // Get ALICE top node
262 nTop=gAlice->GetGeometry()->GetNode("alice");
266 rl = fTPCParam->GetInnerRadiusLow();
267 ru = fTPCParam->GetInnerRadiusUp();
271 sprintf(name,"LS%2.2d",i);
273 sprintf(title,"TPC low sector %3d",i);
276 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
277 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
278 tubs->SetNumberOfDivisions(1);
280 nNode = new TNode(name,title,name,0,0,0,"");
281 nNode->SetLineColor(kColorTPC);
287 rl = fTPCParam->GetOuterRadiusLow();
288 ru = fTPCParam->GetOuterRadiusUp();
291 sprintf(name,"US%2.2d",i);
293 sprintf(title,"TPC upper sector %d",i);
295 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
296 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
297 tubs->SetNumberOfDivisions(1);
299 nNode = new TNode(name,title,name,0,0,0,"");
300 nNode->SetLineColor(kColorTPC);
306 //_____________________________________________________________________________
307 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
310 // Calculate distance from TPC to mouse on the display
316 void AliTPC::Clusters2Tracks()
318 //-----------------------------------------------------------------
319 // This is a track finder.
320 //-----------------------------------------------------------------
321 Error("Clusters2Tracks",
322 "Dummy function ! Call AliTPCtracker::Clusters2Tracks(...) instead !");
325 //_____________________________________________________________________________
326 void AliTPC::CreateMaterials()
328 //-----------------------------------------------
329 // Create Materials for for TPC simulations
330 //-----------------------------------------------
332 //-----------------------------------------------------------------
333 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
334 //-----------------------------------------------------------------
336 Int_t iSXFLD=gAlice->Field()->Integ();
337 Float_t sXMGMX=gAlice->Field()->Max();
339 Float_t amat[5]; // atomic numbers
340 Float_t zmat[5]; // z
341 Float_t wmat[5]; // proportions
347 //***************** Gases *************************
349 //-------------------------------------------------
351 //-------------------------------------------------
362 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
372 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
375 //--------------------------------------------------------------
377 //--------------------------------------------------------------
394 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
396 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
411 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
413 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
429 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
431 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
433 //----------------------------------------------------------------
434 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
435 //----------------------------------------------------------------
441 Float_t rho,absl,X0,buf[1];
445 for(nc = 0;nc<fNoComp;nc++)
448 // retrive material constants
450 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
455 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
457 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
458 density += fMixtProp[nc]*rho; // density of the mixture
462 // mixture proportions by weight!
464 for(nc = 0;nc<fNoComp;nc++)
467 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
469 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
470 apure[nnc] : amol[nnc])/am;
474 // Drift gases 1 - nonsensitive, 2 - sensitive
476 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
477 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
486 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
489 //----------------------------------------------------------------------
491 //----------------------------------------------------------------------
513 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
535 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
553 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
571 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
573 // SiO2 - used later for the glass fiber
585 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
594 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
603 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
612 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
630 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
649 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
651 // Epoxy - C14 H20 O3
668 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
676 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
680 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,X0,absl,buf,nbuf);
684 wmat[0]=0.644; // by weight!
687 density=0.5*(1.25+2.265);
689 AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
693 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
695 wmat[0]=0.725; // by weight!
700 AliMixture(39,"G10",amat,zmat,density,2,wmat);
705 //----------------------------------------------------------
706 // tracking media for gases
707 //----------------------------------------------------------
709 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
710 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
711 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
712 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
714 //-----------------------------------------------------------
715 // tracking media for solids
716 //-----------------------------------------------------------
718 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
719 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
720 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
721 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
722 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
723 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
724 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
725 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
726 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
727 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
728 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
729 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
733 void AliTPC::GenerNoise(Int_t tablesize)
736 //Generate table with noise
743 if (fNoiseTable) delete[] fNoiseTable;
744 fNoiseTable = new Float_t[tablesize];
745 fNoiseDepth = tablesize;
746 fCurrentNoise =0; //!index of the noise in the noise table
748 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
749 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
752 Float_t AliTPC::GetNoise()
754 // get noise from table
755 // if ((fCurrentNoise%10)==0)
756 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
757 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
758 return fNoiseTable[fCurrentNoise++];
759 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
763 Bool_t AliTPC::IsSectorActive(Int_t sec)
766 // check if the sector is active
767 if (!fActiveSectors) return kTRUE;
768 else return fActiveSectors[sec];
771 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
773 // activate interesting sectors
774 SetTreeAddress();//just for security
775 if (fActiveSectors) delete [] fActiveSectors;
776 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
777 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
778 for (Int_t i=0;i<n;i++)
779 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
783 void AliTPC::SetActiveSectors(Int_t flag)
786 // activate sectors which were hitted by tracks
788 SetTreeAddress();//just for security
789 if (fHitType==0) return; // if Clones hit - not short volume ID information
790 if (fActiveSectors) delete [] fActiveSectors;
791 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
793 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
796 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
800 Fatal("SetActiveSectors","Can not find TreeH in folder");
803 if (fHitType>1) branch = TreeH()->GetBranch("TPC2");
804 else branch = TreeH()->GetBranch("TPC");
805 Stat_t ntracks = TreeH()->GetEntries();
806 // loop over all hits
807 cout<<"\nAliTPC::SetActiveSectors(): Got "<<ntracks<<" tracks\n";
809 for(Int_t track=0;track<ntracks;track++)
813 if (fTrackHits && fHitType&4) {
814 TBranch * br1 = TreeH()->GetBranch("fVolumes");
815 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
816 br1->GetEvent(track);
817 br2->GetEvent(track);
818 Int_t *volumes = fTrackHits->GetVolumes();
819 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
820 fActiveSectors[volumes[j]]=kTRUE;
824 if (fTrackHitsOld && fHitType&2) {
825 TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
827 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
828 for (UInt_t j=0;j<ar->GetSize();j++){
829 fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
839 void AliTPC::Digits2Clusters(Int_t /*eventnumber*/)
841 //-----------------------------------------------------------------
842 // This is a simple cluster finder.
843 //-----------------------------------------------------------------
844 Error("Digits2Clusters",
845 "Dummy function ! Call AliTPCclusterer::Digits2Clusters(...) instead !");
848 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
849 extern Double_t SigmaZ2(Double_t, Double_t);
850 //_____________________________________________________________________________
851 void AliTPC::Hits2Clusters(Int_t /*eventn*/)
853 //--------------------------------------------------------
854 // TPC simple cluster generator from hits
855 // obtained from the TPC Fast Simulator
856 // The point errors are taken from the parametrization
857 //--------------------------------------------------------
859 //-----------------------------------------------------------------
860 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
861 //-----------------------------------------------------------------
862 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
863 // Jouri.Belikov@cern.ch
864 //----------------------------------------------------------------
866 /////////////////////////////////////////////////////////////////////////////
868 //---------------------------------------------------------------------
869 // ALICE TPC Cluster Parameters
870 //--------------------------------------------------------------------
874 // Cluster width in rphi
875 const Float_t kACrphi=0.18322;
876 const Float_t kBCrphi=0.59551e-3;
877 const Float_t kCCrphi=0.60952e-1;
878 // Cluster width in z
879 const Float_t kACz=0.19081;
880 const Float_t kBCz=0.55938e-3;
881 const Float_t kCCz=0.30428;
885 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
889 //if(fDefaults == 0) SetDefaults();
891 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
893 TParticle *particle; // pointer to a given particle
894 AliTPChit *tpcHit; // pointer to a sigle TPC hit
898 Float_t pl,pt,tanth,rpad,ratio;
901 //---------------------------------------------------------------
902 // Get the access to the tracks
903 //---------------------------------------------------------------
908 Fatal("Hits2Clusters","Can not find TreeH in folder");
913 Stat_t ntracks = tH->GetEntries();
915 //Switch to the output file
917 if (fLoader->TreeR() == 0x0) fLoader->MakeTree("R");
919 cout<<"fTPCParam->GetTitle() = "<<fTPCParam->GetTitle()<<endl;
921 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
923 //fTPCParam->Write(fTPCParam->GetTitle());
925 AliTPCClustersArray carray;
926 carray.Setup(fTPCParam);
927 carray.SetClusterType("AliTPCcluster");
928 carray.MakeTree(fLoader->TreeR());
930 Int_t nclusters=0; //cluster counter
932 //------------------------------------------------------------
933 // Loop over all sectors (72 sectors for 20 deg
934 // segmentation for both lower and upper sectors)
935 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
936 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
938 // First cluster for sector 0 starts at "0"
939 //------------------------------------------------------------
941 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
943 fTPCParam->AdjustCosSin(isec,cph,sph);
945 //------------------------------------------------------------
947 //------------------------------------------------------------
949 for(Int_t track=0;track<ntracks;track++){
953 // Get number of the TPC hits
955 tpcHit = (AliTPChit*)FirstHit(-1);
961 if (tpcHit->fQ == 0.) {
962 tpcHit = (AliTPChit*) NextHit();
963 continue; //information about track (I.Belikov)
965 sector=tpcHit->fSector; // sector number
968 tpcHit = (AliTPChit*) NextHit();
971 ipart=tpcHit->Track();
972 particle=gAlice->GetMCApp()->Particle(ipart);
975 if(pt < 1.e-9) pt=1.e-9;
977 tanth = TMath::Abs(tanth);
978 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
979 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
981 // space-point resolutions
983 sigmaRphi=SigmaY2(rpad,tanth,pt);
984 sigmaZ =SigmaZ2(rpad,tanth );
988 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
989 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
991 // temporary protection
993 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
994 if(sigmaZ < 0.) sigmaZ=0.4e-3;
995 if(clRphi < 0.) clRphi=2.5e-3;
996 if(clZ < 0.) clZ=2.5e-5;
1001 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
1002 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1004 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
1005 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
1006 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
1007 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
1008 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1009 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
1010 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
1011 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
1012 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
1013 xyz[2]=sigmaRphi; // fSigmaY2
1014 xyz[3]=sigmaZ; // fSigmaZ2
1015 xyz[4]=tpcHit->fQ; // q
1017 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
1018 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
1020 Int_t tracks[3]={tpcHit->Track(), -1, -1};
1021 AliTPCcluster cluster(tracks,xyz);
1023 clrow->InsertCluster(&cluster); nclusters++;
1025 tpcHit = (AliTPChit*)NextHit();
1028 } // end of loop over hits
1030 } // end of loop over tracks
1032 Int_t nrows=fTPCParam->GetNRow(isec);
1033 for (Int_t irow=0; irow<nrows; irow++) {
1034 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
1035 if (!clrow) continue;
1036 carray.StoreRow(isec,irow);
1037 carray.ClearRow(isec,irow);
1040 } // end of loop over sectors
1042 cerr<<"Number of made clusters : "<<nclusters<<" \n";
1043 fLoader->WriteRecPoints("OVERWRITE");
1046 } // end of function
1048 //_________________________________________________________________
1049 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1051 //--------------------------------------------------------
1052 //calculate exact cross point of track and given pad row
1053 //resulting values are expressed in "digit" coordinata
1054 //--------------------------------------------------------
1056 //-----------------------------------------------------------------
1057 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1058 //-----------------------------------------------------------------
1060 if (fClustersArray==0){
1064 TParticle *particle; // pointer to a given particle
1065 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1066 // Int_t sector,nhits;
1068 const Int_t kcmaxhits=30000;
1069 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
1070 AliTPCFastVector & xxx = *xxxx;
1071 Int_t maxhits = kcmaxhits;
1072 //construct array for each padrow
1073 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
1074 fClustersArray->CreateRow(isec,i);
1076 //---------------------------------------------------------------
1077 // Get the access to the tracks
1078 //---------------------------------------------------------------
1080 TTree *tH = TreeH();
1083 Fatal("Hits2Clusters","Can not find TreeH in folder");
1088 Stat_t ntracks = tH->GetEntries();
1089 Int_t npart = gAlice->GetMCApp()->GetNtrack();
1092 if (fHitType>1) branch = tH->GetBranch("TPC2");
1093 else branch = tH->GetBranch("TPC");
1095 //------------------------------------------------------------
1097 //------------------------------------------------------------
1099 for(Int_t track=0;track<ntracks;track++){
1100 Bool_t isInSector=kTRUE;
1102 isInSector = TrackInVolume(isec,track);
1103 if (!isInSector) continue;
1105 branch->GetEntry(track); // get next track
1107 // Get number of the TPC hits and a pointer
1111 Int_t currentIndex=0;
1112 Int_t lastrow=-1; //last writen row
1116 tpcHit = (AliTPChit*)FirstHit(-1);
1119 Int_t sector=tpcHit->fSector; // sector number
1121 tpcHit = (AliTPChit*) NextHit();
1125 ipart=tpcHit->Track();
1126 if (ipart<npart) particle=gAlice->GetMCApp()->Particle(ipart);
1130 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1131 Int_t index[3]={1,isec,0};
1132 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1133 if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
1134 if (lastrow<0) lastrow=currentrow;
1135 if (currentrow==lastrow){
1136 if ( currentIndex>=maxhits){
1138 xxx.ResizeTo(4*maxhits);
1140 xxx(currentIndex*4)=x[0];
1141 xxx(currentIndex*4+1)=x[1];
1142 xxx(currentIndex*4+2)=x[2];
1143 xxx(currentIndex*4+3)=tpcHit->fQ;
1147 if (currentIndex>2){
1159 for (Int_t index=0;index<currentIndex;index++){
1161 x=x2=x3=x4=xxx(index*4);
1169 sumy+=xxx(index*4+1);
1170 sumxy+=xxx(index*4+1)*x;
1171 sumx2y+=xxx(index*4+1)*x2;
1172 sumz+=xxx(index*4+2);
1173 sumxz+=xxx(index*4+2)*x;
1174 sumx2z+=xxx(index*4+2)*x2;
1175 sumq+=xxx(index*4+3);
1177 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1178 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1179 sumx2*(sumx*sumx3-sumx2*sumx2);
1181 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1182 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1183 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1184 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1186 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1187 sumx2*(sumx*sumx2y-sumx2*sumxy);
1188 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1189 sumx2*(sumx*sumx2z-sumx2*sumxz);
1191 if (TMath::Abs(det)<0.00001){
1192 tpcHit = (AliTPChit*)NextHit();
1196 Float_t y=detay/det+centralPad;
1197 Float_t z=detaz/det;
1198 Float_t by=detby/det; //y angle
1199 Float_t bz=detbz/det; //z angle
1200 sumy/=Float_t(currentIndex);
1201 sumz/=Float_t(currentIndex);
1203 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1205 AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
1206 // AliComplexCluster cl;
1212 cl->fTracks[0]=ipart;
1216 } //end of calculating cluster for given row
1219 tpcHit = (AliTPChit*)NextHit();
1220 } // end of loop over hits
1221 } // end of loop over tracks
1222 //write padrows to tree
1223 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1224 fClustersArray->StoreRow(isec,ii);
1225 fClustersArray->ClearRow(isec,ii);
1233 //______________________________________________________________________
1234 AliDigitizer* AliTPC::CreateDigitizer(AliRunDigitizer* manager)
1236 return new AliTPCDigitizer(manager);
1239 void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
1241 //create digits from summable digits
1242 GenerNoise(500000); //create teble with noise
1244 //conect tree with sSDigits
1245 TTree *t = fLoader->TreeS();
1249 fLoader->LoadSDigits("READ");
1250 t = fLoader->TreeS();
1253 Error("SDigits2Digits2","Can not get input TreeS");
1258 if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
1260 AliSimDigits digarr, *dummy=&digarr;
1261 TBranch* sdb = t->GetBranch("Segment");
1264 Error("SDigits2Digits2","Can not find branch with segments in TreeS.");
1268 sdb->SetAddress(&dummy);
1270 Stat_t nentries = t->GetEntries();
1272 // set zero suppression
1274 fTPCParam->SetZeroSup(2);
1276 // get zero suppression
1278 Int_t zerosup = fTPCParam->GetZeroSup();
1280 //make tree with digits
1282 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1283 arr->SetClass("AliSimDigits");
1284 arr->Setup(fTPCParam);
1285 arr->MakeTree(fLoader->TreeD());
1287 AliTPCParam * par = fTPCParam;
1289 //Loop over segments of the TPC
1291 for (Int_t n=0; n<nentries; n++) {
1294 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1295 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
1298 if (!IsSectorActive(sec))
1300 // cout<<n<<" NOT Active \n";
1305 // cout<<n<<" Active \n";
1307 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1308 Int_t nrows = digrow->GetNRows();
1309 Int_t ncols = digrow->GetNCols();
1311 digrow->ExpandBuffer();
1312 digarr.ExpandBuffer();
1313 digrow->ExpandTrackBuffer();
1314 digarr.ExpandTrackBuffer();
1317 Short_t * pamp0 = digarr.GetDigits();
1318 Int_t * ptracks0 = digarr.GetTracks();
1319 Short_t * pamp1 = digrow->GetDigits();
1320 Int_t * ptracks1 = digrow->GetTracks();
1321 Int_t nelems =nrows*ncols;
1322 Int_t saturation = fTPCParam->GetADCSat();
1323 //use internal structure of the AliDigits - for speed reason
1324 //if you cahnge implementation
1325 //of the Alidigits - it must be rewriten -
1326 for (Int_t i= 0; i<nelems; i++){
1327 // Float_t q = *pamp0;
1328 //q/=16.; //conversion faktor
1329 //Float_t noise= GetNoise();
1331 //q= TMath::Nint(q);
1332 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1334 if (q>saturation) q=saturation;
1336 //if (ptracks0[0]==0)
1339 ptracks1[0]=ptracks0[0];
1340 ptracks1[nelems]=ptracks0[nelems];
1341 ptracks1[2*nelems]=ptracks0[2*nelems];
1349 arr->StoreRow(sec,row);
1350 arr->ClearRow(sec,row);
1351 // cerr<<sec<<"\t"<<row<<"\n";
1356 fLoader->WriteDigits("OVERWRITE");
1360 //__________________________________________________________________
1361 void AliTPC::SetDefaults(){
1364 cerr<<"Setting default parameters...\n";
1366 // Set response functions
1369 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1371 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1373 printf("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...\n");
1375 param = new AliTPCParamSR();
1378 param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1381 printf("No TPC parameters found\n");
1386 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1387 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1388 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1389 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1390 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1391 rf->SetOffset(3*param->GetZSigma());
1394 TDirectory *savedir=gDirectory;
1395 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1397 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1402 prfinner->Read("prf_07504_Gati_056068_d02");
1403 //PH Set different names
1404 s=prfinner->GetGRF()->GetName();
1406 prfinner->GetGRF()->SetName(s.Data());
1408 prfouter1->Read("prf_10006_Gati_047051_d03");
1409 s=prfouter1->GetGRF()->GetName();
1411 prfouter1->GetGRF()->SetName(s.Data());
1413 prfouter2->Read("prf_15006_Gati_047051_d03");
1414 s=prfouter2->GetGRF()->GetName();
1416 prfouter2->GetGRF()->SetName(s.Data());
1421 param->SetInnerPRF(prfinner);
1422 param->SetOuter1PRF(prfouter1);
1423 param->SetOuter2PRF(prfouter2);
1424 param->SetTimeRF(rf);
1434 //__________________________________________________________________
1435 void AliTPC::Hits2Digits()
1437 fLoader->LoadHits("read");
1438 fLoader->LoadDigits("recreate");
1439 AliRunLoader* runLoader = fLoader->GetRunLoader();
1441 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1442 runLoader->GetEvent(iEvent);
1444 Hits2Digits(iEvent);
1447 fLoader->UnloadHits();
1448 fLoader->UnloadDigits();
1450 //__________________________________________________________________
1451 void AliTPC::Hits2Digits(Int_t eventnumber)
1453 //----------------------------------------------------
1454 // Loop over all sectors for a single event
1455 //----------------------------------------------------
1456 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1457 rl->GetEvent(eventnumber);
1458 if (fLoader->TreeH() == 0x0)
1460 if(fLoader->LoadHits())
1462 Error("Hits2Digits","Can not load hits.");
1467 if (fLoader->TreeD() == 0x0 )
1469 fLoader->MakeTree("D");
1470 if (fLoader->TreeD() == 0x0 )
1472 Error("Hits2Digits","Can not get TreeD");
1477 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1478 GenerNoise(500000); //create teble with noise
1480 //setup TPCDigitsArray
1482 if(GetDigitsArray()) delete GetDigitsArray();
1484 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1485 arr->SetClass("AliSimDigits");
1486 arr->Setup(fTPCParam);
1488 arr->MakeTree(fLoader->TreeD());
1489 SetDigitsArray(arr);
1491 fDigitsSwitch=0; // standard digits
1493 cerr<<"Digitizing TPC -- normal digits...\n";
1495 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1496 if (IsSectorActive(isec))
1498 if (fDebug) Info("Hits2Digits","Sector %d is active.",isec);
1499 Hits2DigitsSector(isec);
1503 if (fDebug) Info("Hits2Digits","Sector %d is NOT active.",isec);
1506 fLoader->WriteDigits("OVERWRITE");
1508 //this line prevents the crash in the similar one
1509 //on the beginning of this method
1510 //destructor attempts to reset the tree, which is deleted by the loader
1511 //need to be redesign
1512 if(GetDigitsArray()) delete GetDigitsArray();
1513 SetDigitsArray(0x0);
1517 //__________________________________________________________________
1518 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1521 //-----------------------------------------------------------
1522 // summable digits - 16 bit "ADC", no noise, no saturation
1523 //-----------------------------------------------------------
1525 //----------------------------------------------------
1526 // Loop over all sectors for a single event
1527 //----------------------------------------------------
1528 // AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1530 AliRunLoader* rl = fLoader->GetRunLoader();
1532 rl->GetEvent(eventnumber);
1533 if (fLoader->TreeH() == 0x0)
1535 if(fLoader->LoadHits())
1537 Error("Hits2Digits","Can not load hits.");
1544 if (fLoader->TreeS() == 0x0 )
1546 fLoader->MakeTree("S");
1549 if(fDefaults == 0) SetDefaults();
1551 GenerNoise(500000); //create table with noise
1552 //setup TPCDigitsArray
1554 if(GetDigitsArray()) delete GetDigitsArray();
1557 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1558 arr->SetClass("AliSimDigits");
1559 arr->Setup(fTPCParam);
1560 arr->MakeTree(fLoader->TreeS());
1562 SetDigitsArray(arr);
1564 cerr<<"Digitizing TPC -- summable digits...\n";
1566 fDigitsSwitch=1; // summable digits
1568 // set zero suppression to "0"
1570 fTPCParam->SetZeroSup(0);
1572 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1573 if (IsSectorActive(isec))
1575 // cout<<"Sector "<<isec<<" is active\n";
1576 Hits2DigitsSector(isec);
1579 fLoader->WriteSDigits("OVERWRITE");
1581 //this line prevents the crash in the similar one
1582 //on the beginning of this method
1583 //destructor attempts to reset the tree, which is deleted by the loader
1584 //need to be redesign
1585 if(GetDigitsArray()) delete GetDigitsArray();
1586 SetDigitsArray(0x0);
1588 //__________________________________________________________________
1590 void AliTPC::Hits2SDigits()
1593 //-----------------------------------------------------------
1594 // summable digits - 16 bit "ADC", no noise, no saturation
1595 //-----------------------------------------------------------
1597 fLoader->LoadHits("read");
1598 fLoader->LoadSDigits("recreate");
1599 AliRunLoader* runLoader = fLoader->GetRunLoader();
1601 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1602 runLoader->GetEvent(iEvent);
1605 Hits2SDigits2(iEvent);
1608 fLoader->UnloadHits();
1609 fLoader->UnloadSDigits();
1611 //_____________________________________________________________________________
1613 void AliTPC::Hits2DigitsSector(Int_t isec)
1615 //-------------------------------------------------------------------
1616 // TPC conversion from hits to digits.
1617 //-------------------------------------------------------------------
1619 //-----------------------------------------------------------------
1620 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1621 //-----------------------------------------------------------------
1623 //-------------------------------------------------------
1624 // Get the access to the track hits
1625 //-------------------------------------------------------
1627 // check if the parameters are set - important if one calls this method
1628 // directly, not from the Hits2Digits
1630 if(fDefaults == 0) SetDefaults();
1632 TTree *tH = TreeH(); // pointer to the hits tree
1635 Fatal("Hits2DigitsSector","Can not find TreeH in folder");
1639 Stat_t ntracks = tH->GetEntries();
1643 //-------------------------------------------
1644 // Only if there are any tracks...
1645 //-------------------------------------------
1649 //printf("*** Processing sector number %d ***\n",isec);
1651 Int_t nrows =fTPCParam->GetNRow(isec);
1653 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1655 MakeSector(isec,nrows,tH,ntracks,row);
1657 //--------------------------------------------------------
1658 // Digitize this sector, row by row
1659 // row[i] is the pointer to the TObjArray of AliTPCFastVectors,
1660 // each one containing electrons accepted on this
1661 // row, assigned into tracks
1662 //--------------------------------------------------------
1666 if (fDigitsArray->GetTree()==0)
1668 Fatal("Hits2DigitsSector","Tree not set in fDigitsArray");
1671 for (i=0;i<nrows;i++){
1673 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1675 DigitizeRow(i,isec,row);
1677 fDigitsArray->StoreRow(isec,i);
1679 Int_t ndig = dig->GetDigitSize();
1682 printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1684 fDigitsArray->ClearRow(isec,i);
1687 } // end of the sector digitization
1689 for(i=0;i<nrows+2;i++){
1694 delete [] row; // delete the array of pointers to TObjArray-s
1698 } // end of Hits2DigitsSector
1701 //_____________________________________________________________________________
1702 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1704 //-----------------------------------------------------------
1705 // Single row digitization, coupling from the neighbouring
1706 // rows taken into account
1707 //-----------------------------------------------------------
1709 //-----------------------------------------------------------------
1710 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1711 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1712 //-----------------------------------------------------------------
1715 Float_t zerosup = fTPCParam->GetZeroSup();
1716 // Int_t nrows =fTPCParam->GetNRow(isec);
1717 fCurrentIndex[1]= isec;
1720 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1721 Int_t nofTbins = fTPCParam->GetMaxTBin();
1722 Int_t indexRange[4];
1724 // Integrated signal for this row
1725 // and a single track signal
1728 AliTPCFastMatrix *m1 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // integrated
1729 AliTPCFastMatrix *m2 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // single
1731 AliTPCFastMatrix &total = *m1;
1733 // Array of pointers to the label-signal list
1735 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1736 Float_t **pList = new Float_t* [nofDigits];
1740 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1744 //Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1745 //Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1748 for (Int_t row= row1;row<=row2;row++){
1749 Int_t nTracks= rows[row]->GetEntries();
1750 for (i1=0;i1<nTracks;i1++){
1751 fCurrentIndex[2]= row;
1752 fCurrentIndex[3]=irow+1;
1754 m2->Zero(); // clear single track signal matrix
1755 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1756 GetList(trackLabel,nofPads,m2,indexRange,pList);
1758 else GetSignal(rows[row],i1,0,m1,indexRange);
1764 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1766 Float_t fzerosup = zerosup+0.5;
1767 for(Int_t it=0;it<nofTbins;it++){
1768 Float_t *pq = &(total.UncheckedAt(0,it));
1769 for(Int_t ip=0;ip<nofPads;ip++){
1773 if(fDigitsSwitch == 0){
1775 if(q <=fzerosup) continue; // do not fill zeros
1777 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1782 if(q <= 0.) continue; // do not fill zeros
1783 if(q>2000.) q=2000.;
1789 // "real" signal or electronic noise (list = -1)?
1792 for(Int_t j1=0;j1<3;j1++){
1793 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1798 <A NAME="AliDigits"></A>
1799 using of AliDigits object
1802 dig->SetDigitFast((Short_t)q,it,ip);
1803 if (fDigitsArray->IsSimulated())
1805 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1806 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1807 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1811 } // end of loop over time buckets
1812 } // end of lop over pads
1815 // This row has been digitized, delete nonused stuff
1818 for(lp=0;lp<nofDigits;lp++){
1819 if(pList[lp]) delete [] pList[lp];
1828 } // end of DigitizeRow
1830 //_____________________________________________________________________________
1832 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1833 AliTPCFastMatrix *m1, AliTPCFastMatrix *m2,Int_t *indexRange)
1836 //---------------------------------------------------------------
1837 // Calculates 2-D signal (pad,time) for a single track,
1838 // returns a pointer to the signal matrix and the track label
1839 // No digitization is performed at this level!!!
1840 //---------------------------------------------------------------
1842 //-----------------------------------------------------------------
1843 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1844 // Modified: Marian Ivanov
1845 //-----------------------------------------------------------------
1847 AliTPCFastVector *tv;
1849 tv = (AliTPCFastVector*)p1->At(ntr); // pointer to a track
1850 AliTPCFastVector &v = *tv;
1852 Float_t label = v(0);
1853 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1)-1)/2;
1855 Int_t nElectrons = (tv->GetNrows()-1)/4;
1856 indexRange[0]=9999; // min pad
1857 indexRange[1]=-1; // max pad
1858 indexRange[2]=9999; //min time
1859 indexRange[3]=-1; // max time
1861 AliTPCFastMatrix &signal = *m1;
1862 AliTPCFastMatrix &total = *m2;
1864 // Loop over all electrons
1866 for(Int_t nel=0; nel<nElectrons; nel++){
1868 Float_t aval = v(idx+4);
1869 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1870 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1871 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
1873 Int_t *index = fTPCParam->GetResBin(0);
1874 Float_t *weight = & (fTPCParam->GetResWeight(0));
1876 if (n>0) for (Int_t i =0; i<n; i++){
1877 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1880 Int_t time=index[2];
1881 Float_t qweight = *(weight)*eltoadcfac;
1883 if (m1!=0) signal.UncheckedAt(pad,time)+=qweight;
1884 total.UncheckedAt(pad,time)+=qweight;
1885 if (indexRange[0]>pad) indexRange[0]=pad;
1886 if (indexRange[1]<pad) indexRange[1]=pad;
1887 if (indexRange[2]>time) indexRange[2]=time;
1888 if (indexRange[3]<time) indexRange[3]=time;
1895 } // end of loop over electrons
1897 return label; // returns track label when finished
1900 //_____________________________________________________________________________
1901 void AliTPC::GetList(Float_t label,Int_t np,AliTPCFastMatrix *m,
1902 Int_t *indexRange, Float_t **pList)
1904 //----------------------------------------------------------------------
1905 // Updates the list of tracks contributing to digits for a given row
1906 //----------------------------------------------------------------------
1908 //-----------------------------------------------------------------
1909 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1910 //-----------------------------------------------------------------
1912 AliTPCFastMatrix &signal = *m;
1914 // lop over nonzero digits
1916 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1917 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1920 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1922 if(signal(ip,it)<0.5) continue;
1925 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1927 if(!pList[globalIndex]){
1930 // Create new list (6 elements - 3 signals and 3 labels),
1933 pList[globalIndex] = new Float_t [6];
1937 *pList[globalIndex] = -1.;
1938 *(pList[globalIndex]+1) = -1.;
1939 *(pList[globalIndex]+2) = -1.;
1940 *(pList[globalIndex]+3) = -1.;
1941 *(pList[globalIndex]+4) = -1.;
1942 *(pList[globalIndex]+5) = -1.;
1945 *pList[globalIndex] = label;
1946 *(pList[globalIndex]+3) = signal(ip,it);
1950 // check the signal magnitude
1952 Float_t highest = *(pList[globalIndex]+3);
1953 Float_t middle = *(pList[globalIndex]+4);
1954 Float_t lowest = *(pList[globalIndex]+5);
1957 // compare the new signal with already existing list
1960 if(signal(ip,it)<lowest) continue; // neglect this track
1964 if (signal(ip,it)>highest){
1965 *(pList[globalIndex]+5) = middle;
1966 *(pList[globalIndex]+4) = highest;
1967 *(pList[globalIndex]+3) = signal(ip,it);
1969 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1970 *(pList[globalIndex]+1) = *pList[globalIndex];
1971 *pList[globalIndex] = label;
1973 else if (signal(ip,it)>middle){
1974 *(pList[globalIndex]+5) = middle;
1975 *(pList[globalIndex]+4) = signal(ip,it);
1977 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1978 *(pList[globalIndex]+1) = label;
1981 *(pList[globalIndex]+5) = signal(ip,it);
1982 *(pList[globalIndex]+2) = label;
1986 } // end of loop over pads
1987 } // end of loop over time bins
1992 //___________________________________________________________________
1993 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1994 Stat_t ntracks,TObjArray **row)
1997 //-----------------------------------------------------------------
1998 // Prepares the sector digitization, creates the vectors of
1999 // tracks for each row of this sector. The track vector
2000 // contains the track label and the position of electrons.
2001 //-----------------------------------------------------------------
2003 //-----------------------------------------------------------------
2004 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2005 //-----------------------------------------------------------------
2007 Float_t gasgain = fTPCParam->GetGasGain();
2011 AliTPChit *tpcHit; // pointer to a sigle TPC hit
2014 if (fHitType>1) branch = TH->GetBranch("TPC2");
2015 else branch = TH->GetBranch("TPC");
2018 //----------------------------------------------
2019 // Create TObjArray-s, one for each row,
2020 // each TObjArray will store the AliTPCFastVectors
2021 // of electrons, one AliTPCFastVectors per each track.
2022 //----------------------------------------------
2024 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
2025 AliTPCFastVector **tracks = new AliTPCFastVector* [nrows+2]; //pointers to the track vectors
2027 for(i=0; i<nrows+2; i++){
2028 row[i] = new TObjArray;
2035 //--------------------------------------------------------------------
2036 // Loop over tracks, the "track" contains the full history
2037 //--------------------------------------------------------------------
2039 Int_t previousTrack,currentTrack;
2040 previousTrack = -1; // nothing to store so far!
2042 for(Int_t track=0;track<ntracks;track++){
2043 Bool_t isInSector=kTRUE;
2045 isInSector = TrackInVolume(isec,track);
2046 if (!isInSector) continue;
2048 branch->GetEntry(track); // get next track
2052 tpcHit = (AliTPChit*)FirstHit(-1);
2054 //--------------------------------------------------------------
2056 //--------------------------------------------------------------
2061 Int_t sector=tpcHit->fSector; // sector number
2063 tpcHit = (AliTPChit*) NextHit();
2067 currentTrack = tpcHit->Track(); // track number
2070 if(currentTrack != previousTrack){
2072 // store already filled fTrack
2074 for(i=0;i<nrows+2;i++){
2075 if(previousTrack != -1){
2076 if(nofElectrons[i]>0){
2077 AliTPCFastVector &v = *tracks[i];
2078 v(0) = previousTrack;
2079 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2080 row[i]->Add(tracks[i]);
2083 delete tracks[i]; // delete empty AliTPCFastVector
2089 tracks[i] = new AliTPCFastVector(481); // AliTPCFastVectors for the next fTrack
2091 } // end of loop over rows
2093 previousTrack=currentTrack; // update track label
2096 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2098 //---------------------------------------------------
2099 // Calculate the electron attachment probability
2100 //---------------------------------------------------
2103 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
2104 /fTPCParam->GetDriftV();
2106 Float_t attProb = fTPCParam->GetAttCoef()*
2107 fTPCParam->GetOxyCont()*time; // fraction!
2109 //-----------------------------------------------
2110 // Loop over electrons
2111 //-----------------------------------------------
2114 for(Int_t nel=0;nel<qI;nel++){
2115 // skip if electron lost due to the attachment
2116 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2121 // protection for the nonphysical avalanche size (10**6 maximum)
2123 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2124 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2127 TransportElectron(xyz,index);
2129 fTPCParam->GetPadRow(xyz,index);
2130 // row 0 - cross talk from the innermost row
2131 // row fNRow+1 cross talk from the outermost row
2132 rowNumber = index[2]+1;
2133 //transform position to local digit coordinates
2134 //relative to nearest pad row
2135 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2137 if (isec <fTPCParam->GetNInnerSector()) {
2138 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2139 y1 = fTPCParam->GetYInner(rowNumber);
2142 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2143 y1 = fTPCParam->GetYOuter(rowNumber);
2145 // gain inefficiency at the wires edges - linear
2148 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));
2150 nofElectrons[rowNumber]++;
2151 //----------------------------------
2152 // Expand vector if necessary
2153 //----------------------------------
2154 if(nofElectrons[rowNumber]>120){
2155 Int_t range = tracks[rowNumber]->GetNrows();
2156 if((nofElectrons[rowNumber])>(range-1)/4){
2158 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
2162 AliTPCFastVector &v = *tracks[rowNumber];
2163 Int_t idx = 4*nofElectrons[rowNumber]-3;
2164 Real_t * position = &(((AliTPCFastVector&)v).UncheckedAt(idx)); //make code faster
2165 memcpy(position,xyz,4*sizeof(Float_t));
2167 } // end of loop over electrons
2169 tpcHit = (AliTPChit*)NextHit();
2171 } // end of loop over hits
2172 } // end of loop over tracks
2175 // store remaining track (the last one) if not empty
2178 for(i=0;i<nrows+2;i++){
2179 if(nofElectrons[i]>0){
2180 AliTPCFastVector &v = *tracks[i];
2181 v(0) = previousTrack;
2182 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2183 row[i]->Add(tracks[i]);
2192 delete [] nofElectrons;
2195 } // end of MakeSector
2198 //_____________________________________________________________________________
2202 // Initialise TPC detector after definition of geometry
2207 printf("\n%s: ",ClassName());
2208 for(i=0;i<35;i++) printf("*");
2209 printf(" TPC_INIT ");
2210 for(i=0;i<35;i++) printf("*");
2211 printf("\n%s: ",ClassName());
2213 for(i=0;i<80;i++) printf("*");
2218 //_____________________________________________________________________________
2219 void AliTPC::MakeBranch(Option_t* option)
2222 // Create Tree branches for the TPC.
2224 if(GetDebug()) Info("MakeBranch","");
2225 Int_t buffersize = 4000;
2226 char branchname[10];
2227 sprintf(branchname,"%s",GetName());
2229 const char *h = strstr(option,"H");
2231 if ( h && (fHitType<=1) && (fHits == 0x0)) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2233 AliDetector::MakeBranch(option);
2235 const char *d = strstr(option,"D");
2237 if (fDigits && fLoader->TreeD() && d)
2239 MakeBranchInTree(gAlice->TreeD(), branchname, &fDigits, buffersize, 0);
2242 if (fHitType>1) MakeBranch2(option,0); // MI change 14.09.2000
2245 //_____________________________________________________________________________
2246 void AliTPC::ResetDigits()
2249 // Reset number of digits and the digits array for this detector
2252 if (fDigits) fDigits->Clear();
2255 //_____________________________________________________________________________
2256 void AliTPC::SetSecAL(Int_t sec)
2258 //---------------------------------------------------
2259 // Activate/deactivate selection for lower sectors
2260 //---------------------------------------------------
2262 //-----------------------------------------------------------------
2263 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2264 //-----------------------------------------------------------------
2268 //_____________________________________________________________________________
2269 void AliTPC::SetSecAU(Int_t sec)
2271 //----------------------------------------------------
2272 // Activate/deactivate selection for upper sectors
2273 //---------------------------------------------------
2275 //-----------------------------------------------------------------
2276 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2277 //-----------------------------------------------------------------
2281 //_____________________________________________________________________________
2282 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2284 //----------------------------------------
2285 // Select active lower sectors
2286 //----------------------------------------
2288 //-----------------------------------------------------------------
2289 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2290 //-----------------------------------------------------------------
2300 //_____________________________________________________________________________
2301 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2302 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2303 Int_t s11 , Int_t s12)
2305 //--------------------------------
2306 // Select active upper sectors
2307 //--------------------------------
2309 //-----------------------------------------------------------------
2310 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2311 //-----------------------------------------------------------------
2327 //_____________________________________________________________________________
2328 void AliTPC::SetSens(Int_t sens)
2331 //-------------------------------------------------------------
2332 // Activates/deactivates the sensitive strips at the center of
2333 // the pad row -- this is for the space-point resolution calculations
2334 //-------------------------------------------------------------
2336 //-----------------------------------------------------------------
2337 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2338 //-----------------------------------------------------------------
2344 void AliTPC::SetSide(Float_t side=0.)
2346 // choice of the TPC side
2351 //____________________________________________________________________________
2352 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2353 Float_t p2,Float_t p3)
2356 // gax mixture definition
2370 //_____________________________________________________________________________
2372 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2375 // electron transport taking into account:
2377 // 2.ExB at the wires
2378 // 3. nonisochronity
2380 // xyz and index must be already transformed to system 1
2383 fTPCParam->Transform1to2(xyz,index);
2386 Float_t driftl=xyz[2];
2387 if(driftl<0.01) driftl=0.01;
2388 driftl=TMath::Sqrt(driftl);
2389 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2390 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2391 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2392 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2393 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2397 if (fTPCParam->GetMWPCReadout()==kTRUE){
2398 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2399 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2401 //add nonisochronity (not implemented yet)
2404 ClassImp(AliTPCdigit)
2406 //_____________________________________________________________________________
2407 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2411 // Creates a TPC digit object
2413 fSector = digits[0];
2414 fPadRow = digits[1];
2417 fSignal = digits[4];
2423 //_____________________________________________________________________________
2424 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2428 // Creates a TPC hit object
2438 //________________________________________________________________________
2439 // Additional code because of the AliTPCTrackHitsV2
2441 void AliTPC::MakeBranch2(Option_t *option,const char */*file*/)
2444 // Create a new branch in the current Root Tree
2445 // The branch of fHits is automatically split
2446 // MI change 14.09.2000
2447 if(GetDebug()) Info("MakeBranch2","");
2448 if (fHitType<2) return;
2449 char branchname[10];
2450 sprintf(branchname,"%s2",GetName());
2452 // Get the pointer to the header
2453 const char *cH = strstr(option,"H");
2455 if (fTrackHits && TreeH() && cH && fHitType&4)
2457 if(GetDebug()) Info("MakeBranch2","Making branch for Type 4 Hits");
2458 TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2461 if (fTrackHitsOld && TreeH() && cH && fHitType&2)
2463 if(GetDebug()) Info("MakeBranch2","Making branch for Type 2 Hits");
2464 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2465 TreeH(),fBufferSize,99);
2466 TreeH()->GetListOfBranches()->Add(branch);
2470 void AliTPC::SetTreeAddress()
2472 //Sets tree address for hits
2475 if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2476 AliDetector::SetTreeAddress();
2478 if (fHitType>1) SetTreeAddress2();
2481 void AliTPC::SetTreeAddress2()
2484 // Set branch address for the TrackHits Tree
2486 if(GetDebug()) Info("SetTreeAddress2","");
2489 char branchname[20];
2490 sprintf(branchname,"%s2",GetName());
2492 // Branch address for hit tree
2493 TTree *treeH = TreeH();
2494 if ((treeH)&&(fHitType&4)) {
2495 branch = treeH->GetBranch(branchname);
2498 branch->SetAddress(&fTrackHits);
2499 if (GetDebug()) Info("SetTreeAddress2","fHitType&4 Setting");
2502 if (GetDebug()) Info("SetTreeAddress2","fHitType&4 Failed (can not find branch)");
2505 if ((treeH)&&(fHitType&2)) {
2506 branch = treeH->GetBranch(branchname);
2509 branch->SetAddress(&fTrackHitsOld);
2510 if (GetDebug()) Info("SetTreeAddress2","fHitType&2 Setting");
2512 else if (GetDebug())
2513 Info("SetTreeAddress2","fHitType&2 Failed (can not find branch)");
2515 //set address to TREETR
2517 TTree *treeTR = TreeTR();
2518 if (treeTR && fTrackReferences) {
2519 branch = treeTR->GetBranch(GetName());
2520 if (branch) branch->SetAddress(&fTrackReferences);
2525 void AliTPC::FinishPrimary()
2527 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2528 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2532 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2535 // add hit to the list
2538 int primary = gAlice->GetMCApp()->GetPrimary(track);
2539 gAlice->GetMCApp()->Particle(primary)->SetBit(kKeepBit);
2543 gAlice->GetMCApp()->FlagTrack(track);
2545 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2546 //if (hit->fTrack!=rtrack)
2547 // cout<<"bad track number\n";
2548 if (fTrackHits && fHitType&4)
2549 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2550 hits[1],hits[2],(Int_t)hits[3]);
2551 if (fTrackHitsOld &&fHitType&2 )
2552 fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2553 hits[1],hits[2],(Int_t)hits[3]);
2557 void AliTPC::ResetHits()
2559 if (fHitType&1) AliDetector::ResetHits();
2560 if (fHitType>1) ResetHits2();
2563 void AliTPC::ResetHits2()
2567 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2568 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2572 AliHit* AliTPC::FirstHit(Int_t track)
2574 if (fHitType>1) return FirstHit2(track);
2575 return AliDetector::FirstHit(track);
2577 AliHit* AliTPC::NextHit()
2579 if (fHitType>1) return NextHit2();
2581 return AliDetector::NextHit();
2584 AliHit* AliTPC::FirstHit2(Int_t track)
2587 // Initialise the hit iterator
2588 // Return the address of the first hit for track
2589 // If track>=0 the track is read from disk
2590 // while if track<0 the first hit of the current
2591 // track is returned
2594 gAlice->ResetHits();
2595 TreeH()->GetEvent(track);
2598 if (fTrackHits && fHitType&4) {
2599 fTrackHits->First();
2600 return fTrackHits->GetHit();
2602 if (fTrackHitsOld && fHitType&2) {
2603 fTrackHitsOld->First();
2604 return fTrackHitsOld->GetHit();
2610 AliHit* AliTPC::NextHit2()
2613 //Return the next hit for the current track
2616 if (fTrackHitsOld && fHitType&2) {
2617 fTrackHitsOld->Next();
2618 return fTrackHitsOld->GetHit();
2622 return fTrackHits->GetHit();
2628 void AliTPC::LoadPoints(Int_t)
2632 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2635 if(fHitType==1) AliDetector::LoadPoints(a);
2636 else LoadPoints2(a);
2643 void AliTPC::RemapTrackHitIDs(Int_t *map)
2645 if (!fTrackHits) return;
2647 if (fTrackHitsOld && fHitType&2){
2648 AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2649 for (UInt_t i=0;i<arr->GetSize();i++){
2650 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2651 info->fTrackID = map[info->fTrackID];
2654 if (fTrackHitsOld && fHitType&4){
2655 TClonesArray * arr = fTrackHits->GetArray();;
2656 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2657 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2658 info->fTrackID = map[info->fTrackID];
2663 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2665 //return bool information - is track in given volume
2666 //load only part of the track information
2667 //return true if current track is in volume
2670 if (fTrackHitsOld && fHitType&2) {
2671 TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
2672 br->GetEvent(track);
2673 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2674 for (UInt_t j=0;j<ar->GetSize();j++){
2675 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2679 if (fTrackHits && fHitType&4) {
2680 TBranch * br1 = TreeH()->GetBranch("fVolumes");
2681 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
2682 br2->GetEvent(track);
2683 br1->GetEvent(track);
2684 Int_t *volumes = fTrackHits->GetVolumes();
2685 Int_t nvolumes = fTrackHits->GetNVolumes();
2686 if (!volumes && nvolumes>0) {
2687 printf("Problematic track\t%d\t%d",track,nvolumes);
2690 for (Int_t j=0;j<nvolumes; j++)
2691 if (volumes[j]==id) return kTRUE;
2695 TBranch * br = TreeH()->GetBranch("fSector");
2696 br->GetEvent(track);
2697 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2698 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2705 //_____________________________________________________________________________
2706 void AliTPC::LoadPoints2(Int_t)
2709 // Store x, y, z of all hits in memory
2711 if (fTrackHits == 0 && fTrackHitsOld==0) return;
2714 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2715 if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2717 if (nhits == 0) return;
2718 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2719 if (fPoints == 0) fPoints = new TObjArray(tracks);
2722 Int_t *ntrk=new Int_t[tracks];
2723 Int_t *limi=new Int_t[tracks];
2724 Float_t **coor=new Float_t*[tracks];
2725 for(Int_t i=0;i<tracks;i++) {
2731 AliPoints *points = 0;
2734 Int_t chunk=nhits/4+1;
2736 // Loop over all the hits and store their position
2738 ahit = FirstHit2(-1);
2740 trk=ahit->GetTrack();
2741 if(ntrk[trk]==limi[trk]) {
2743 // Initialise a new track
2744 fp=new Float_t[3*(limi[trk]+chunk)];
2746 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2747 delete [] coor[trk];
2754 fp[3*ntrk[trk] ] = ahit->X();
2755 fp[3*ntrk[trk]+1] = ahit->Y();
2756 fp[3*ntrk[trk]+2] = ahit->Z();
2764 for(trk=0; trk<tracks; ++trk) {
2766 points = new AliPoints();
2767 points->SetMarkerColor(GetMarkerColor());
2768 points->SetMarkerSize(GetMarkerSize());
2769 points->SetDetector(this);
2770 points->SetParticle(trk);
2771 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2772 fPoints->AddAt(points,trk);
2773 delete [] coor[trk];
2783 //_____________________________________________________________________________
2784 void AliTPC::LoadPoints3(Int_t)
2787 // Store x, y, z of all hits in memory
2788 // - only intersection point with pad row
2789 if (fTrackHits == 0) return;
2791 Int_t nhits = fTrackHits->GetEntriesFast();
2792 if (nhits == 0) return;
2793 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2794 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2795 fPoints->Expand(2*tracks);
2798 Int_t *ntrk=new Int_t[tracks];
2799 Int_t *limi=new Int_t[tracks];
2800 Float_t **coor=new Float_t*[tracks];
2801 for(Int_t i=0;i<tracks;i++) {
2807 AliPoints *points = 0;
2810 Int_t chunk=nhits/4+1;
2812 // Loop over all the hits and store their position
2814 ahit = FirstHit2(-1);
2815 //for (Int_t hit=0;hit<nhits;hit++) {
2819 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2820 trk=ahit->GetTrack();
2821 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2822 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2823 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2824 if (currentrow!=lastrow){
2825 lastrow = currentrow;
2826 //later calculate intersection point
2827 if(ntrk[trk]==limi[trk]) {
2829 // Initialise a new track
2830 fp=new Float_t[3*(limi[trk]+chunk)];
2832 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2833 delete [] coor[trk];
2840 fp[3*ntrk[trk] ] = ahit->X();
2841 fp[3*ntrk[trk]+1] = ahit->Y();
2842 fp[3*ntrk[trk]+2] = ahit->Z();
2849 for(trk=0; trk<tracks; ++trk) {
2851 points = new AliPoints();
2852 points->SetMarkerColor(GetMarkerColor()+1);
2853 points->SetMarkerStyle(5);
2854 points->SetMarkerSize(0.2);
2855 points->SetDetector(this);
2856 points->SetParticle(trk);
2857 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2858 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2859 fPoints->AddAt(points,tracks+trk);
2860 delete [] coor[trk];
2871 void AliTPC::FindTrackHitsIntersection(TClonesArray * /*arr*/)
2875 //fill clones array with intersection of current point with the
2880 const Int_t kcmaxhits=30000;
2881 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
2882 AliTPCFastVector & xxx = *xxxx;
2883 Int_t maxhits = kcmaxhits;
2886 AliTPChit * tpcHit=0;
2887 tpcHit = (AliTPChit*)FirstHit2(-1);
2888 Int_t currentIndex=0;
2889 Int_t lastrow=-1; //last writen row
2892 if (tpcHit==0) continue;
2893 sector=tpcHit->fSector; // sector number
2894 ipart=tpcHit->Track();
2898 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2899 Int_t index[3]={1,sector,0};
2900 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2901 if (currentrow<0) continue;
2902 if (lastrow<0) lastrow=currentrow;
2903 if (currentrow==lastrow){
2904 if ( currentIndex>=maxhits){
2906 xxx.ResizeTo(4*maxhits);
2908 xxx(currentIndex*4)=x[0];
2909 xxx(currentIndex*4+1)=x[1];
2910 xxx(currentIndex*4+2)=x[2];
2911 xxx(currentIndex*4+3)=tpcHit->fQ;
2915 if (currentIndex>2){
2927 for (Int_t index=0;index<currentIndex;index++){
2929 x=x2=x3=x4=xxx(index*4);
2937 sumy+=xxx(index*4+1);
2938 sumxy+=xxx(index*4+1)*x;
2939 sumx2y+=xxx(index*4+1)*x2;
2940 sumz+=xxx(index*4+2);
2941 sumxz+=xxx(index*4+2)*x;
2942 sumx2z+=xxx(index*4+2)*x2;
2943 sumq+=xxx(index*4+3);
2945 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2946 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2947 sumx2*(sumx*sumx3-sumx2*sumx2);
2949 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2950 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2951 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2952 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2954 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2955 sumx2*(sumx*sumx2y-sumx2*sumxy);
2956 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2957 sumx2*(sumx*sumx2z-sumx2*sumxz);
2959 Float_t y=detay/det+centralPad;
2960 Float_t z=detaz/det;
2961 Float_t by=detby/det; //y angle
2962 Float_t bz=detbz/det; //z angle
2963 sumy/=Float_t(currentIndex);
2964 sumz/=Float_t(currentIndex);
2966 AliComplexCluster cl;
2972 cl.fTracks[0]=ipart;
2974 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2975 if (row!=0) row->InsertCluster(&cl);
2978 } //end of calculating cluster for given row
2980 } // end of loop over hits
2984 //_______________________________________________________________________________
2985 void AliTPC::Digits2Reco(Int_t firstevent,Int_t lastevent)
2987 // produces rec points from digits and writes them on the root file
2988 // AliTPCclusters.root
2990 TDirectory *cwd = gDirectory;
2993 AliTPCParamSR *dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60");
2995 printf("You are running 2 pad-length geom hits with 3 pad-length geom digits\n");
2997 dig = new AliTPCParamSR();
3001 dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");
3004 printf("No TPC parameters found\n");
3009 cout<<"AliTPC::Digits2Reco: TPC parameteres have been set"<<endl;
3011 if(!gSystem->Getenv("CONFIG_FILE")){
3012 out=TFile::Open("AliTPCclusters.root","recreate");
3018 tmp1=gSystem->Getenv("CONFIG_FILE_PREFIX");
3019 tmp2=gSystem->Getenv("CONFIG_OUTDIR");
3020 sprintf(tmp3,"%s%s/AliTPCclusters.root",tmp1,tmp2);
3021 out=TFile::Open(tmp3,"recreate");
3025 cout<<"AliTPC::Digits2Reco - determination of rec points begins"<<endl;
3028 for(Int_t iev=firstevent;iev<lastevent+1;iev++){
3030 printf("Processing event %d\n",iev);
3031 Digits2Clusters(iev);
3033 cout<<"AliTPC::Digits2Reco - determination of rec points ended"<<endl;
3042 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
3045 fLoader = new AliTPCLoader(GetName(),topfoldername);
3049 ////////////////////////////////////////////////////////////////////////
3050 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
3052 // load TPC paarmeters from a given file or create new if the object
3053 // is not found there
3054 // 12/05/2003 This method should be moved to the AliTPCLoader
3055 // and one has to decide where to store the TPC parameters
3058 sprintf(paramName,"75x40_100x60_150x60");
3059 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
3061 cout<<"TPC parameters "<<paramName<<" found."<<endl;
3063 cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3065 paramTPC = new AliTPCParamSR;
3069 // the older version of parameters can be accessed with this code.
3070 // In some cases, we have old parameters saved in the file but
3071 // digits were created with new parameters, it can be distinguish
3072 // by the name of TPC TreeD. The code here is just for the case
3073 // we would need to compare with old data, uncomment it if needed.
3075 // char paramName[50];
3076 // sprintf(paramName,"75x40_100x60");
3077 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
3079 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
3081 // sprintf(paramName,"75x40_100x60_150x60");
3082 // paramTPC=(AliTPCParam*)in->Get(paramName);
3084 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
3086 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3088 // paramTPC = new AliTPCParamSR;
3096 //____________________________________________________________________________
3097 Double_t SigmaY2(Double_t r, Double_t tgl, Double_t pt)
3100 // Parametrised error of the cluster reconstruction (pad direction)
3103 const Float_t kArphi=0.41818e-2;
3104 const Float_t kBrphi=0.17460e-4;
3105 const Float_t kCrphi=0.30993e-2;
3106 const Float_t kDrphi=0.41061e-3;
3108 pt=TMath::Abs(pt)*1000.;
3110 tgl=TMath::Abs(tgl);
3111 Double_t s=kArphi - kBrphi*r*tgl + kCrphi*x*x + kDrphi*x;
3112 if (s<0.4e-3) s=0.4e-3;
3113 s*=1.3; //Iouri Belikov
3119 //____________________________________________________________________________
3120 Double_t SigmaZ2(Double_t r, Double_t tgl)
3123 // Parametrised error of the cluster reconstruction (drift direction)
3126 const Float_t kAz=0.39614e-2;
3127 const Float_t kBz=0.22443e-4;
3128 const Float_t kCz=0.51504e-1;
3131 tgl=TMath::Abs(tgl);
3132 Double_t s=kAz - kBz*r*tgl + kCz*tgl*tgl;
3133 if (s<0.4e-3) s=0.4e-3;
3134 s*=1.3; //Iouri Belikov