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>
55 #include "AliArrayBranch.h"
56 #include "AliClusters.h"
57 #include "AliComplexCluster.h"
58 #include "AliDigits.h"
60 #include "AliPoints.h"
62 #include "AliRunLoader.h"
63 #include "AliSimDigits.h"
66 #include "AliTPCClustersArray.h"
67 #include "AliTPCClustersRow.h"
68 #include "AliTPCDigitsArray.h"
69 #include "AliTPCLoader.h"
70 #include "AliTPCPRF2D.h"
71 #include "AliTPCParamSR.h"
72 #include "AliTPCRF1D.h"
73 #include "AliTPCTrackHits.h"
74 #include "AliTPCTrackHitsV2.h"
75 #include "AliTPCcluster.h"
76 #include "AliTrackReference.h"
82 //_____________________________________________________________________________
83 // helper class for fast matrix and vector manipulation - no range checking
84 // origin - Marian Ivanov
86 class AliTPCFastMatrix : public TMatrix {
88 AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb);
89 inline Float_t & UncheckedAt(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
90 inline Float_t UncheckedAtFast(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
93 AliTPCFastMatrix::AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb):
94 TMatrix(row_lwb, row_upb,col_lwb,col_upb)
97 //_____________________________________________________________________________
98 class AliTPCFastVector : public TVector {
100 AliTPCFastVector(Int_t size);
101 virtual ~AliTPCFastVector(){;}
102 inline Float_t & UncheckedAt(Int_t index) const {return fElements[index];} //fast acces
106 AliTPCFastVector::AliTPCFastVector(Int_t size):TVector(size){
109 //_____________________________________________________________________________
113 // Default constructor
124 fHitType = 2; //default CONTAINERS - based on ROOT structure
131 //_____________________________________________________________________________
132 AliTPC::AliTPC(const char *name, const char *title)
133 : AliDetector(name,title)
136 // Standard constructor
140 // Initialise arrays of hits and digits
141 fHits = new TClonesArray("AliTPChit", 176);
142 gAlice->GetMCApp()->AddHitList(fHits);
147 fTrackHits = new AliTPCTrackHitsV2;
148 fTrackHits->SetHitPrecision(0.002);
149 fTrackHits->SetStepPrecision(0.003);
150 fTrackHits->SetMaxDistance(100);
152 fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
153 fTrackHitsOld->SetHitPrecision(0.002);
154 fTrackHitsOld->SetStepPrecision(0.003);
155 fTrackHitsOld->SetMaxDistance(100);
162 // Initialise counters
168 // Initialise color attributes
169 SetMarkerColor(kYellow);
172 // Set TPC parameters
176 if (!strcmp(title,"Default")) {
177 fTPCParam = new AliTPCParamSR;
179 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
185 //_____________________________________________________________________________
196 delete fTrackHits; //MI 15.09.2000
197 delete fTrackHitsOld; //MI 10.12.2001
198 if (fNoiseTable) delete [] fNoiseTable;
202 //_____________________________________________________________________________
203 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
206 // Add a hit to the list
208 // TClonesArray &lhits = *fHits;
209 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
211 TClonesArray &lhits = *fHits;
212 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
215 AddHit2(track,vol,hits);
218 //_____________________________________________________________________________
219 void AliTPC::BuildGeometry()
223 // Build TPC ROOT TNode geometry for the event display
228 const int kColorTPC=19;
229 char name[5], title[25];
230 const Double_t kDegrad=TMath::Pi()/180;
231 const Double_t kRaddeg=180./TMath::Pi();
234 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
235 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
237 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
238 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
240 Int_t nLo = fTPCParam->GetNInnerSector()/2;
241 Int_t nHi = fTPCParam->GetNOuterSector()/2;
243 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
244 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
245 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
246 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
249 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
250 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
256 // Get ALICE top node
259 nTop=gAlice->GetGeometry()->GetNode("alice");
263 rl = fTPCParam->GetInnerRadiusLow();
264 ru = fTPCParam->GetInnerRadiusUp();
268 sprintf(name,"LS%2.2d",i);
270 sprintf(title,"TPC low sector %3d",i);
273 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
274 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
275 tubs->SetNumberOfDivisions(1);
277 nNode = new TNode(name,title,name,0,0,0,"");
278 nNode->SetLineColor(kColorTPC);
284 rl = fTPCParam->GetOuterRadiusLow();
285 ru = fTPCParam->GetOuterRadiusUp();
288 sprintf(name,"US%2.2d",i);
290 sprintf(title,"TPC upper sector %d",i);
292 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
293 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
294 tubs->SetNumberOfDivisions(1);
296 nNode = new TNode(name,title,name,0,0,0,"");
297 nNode->SetLineColor(kColorTPC);
303 //_____________________________________________________________________________
304 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
307 // Calculate distance from TPC to mouse on the display
313 void AliTPC::Clusters2Tracks()
315 //-----------------------------------------------------------------
316 // This is a track finder.
317 //-----------------------------------------------------------------
318 Error("Clusters2Tracks",
319 "Dummy function ! Call AliTPCtracker::Clusters2Tracks(...) instead !");
322 //_____________________________________________________________________________
323 void AliTPC::CreateMaterials()
325 //-----------------------------------------------
326 // Create Materials for for TPC simulations
327 //-----------------------------------------------
329 //-----------------------------------------------------------------
330 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
331 //-----------------------------------------------------------------
333 Int_t iSXFLD=gAlice->Field()->Integ();
334 Float_t sXMGMX=gAlice->Field()->Max();
336 Float_t amat[5]; // atomic numbers
337 Float_t zmat[5]; // z
338 Float_t wmat[5]; // proportions
344 //***************** Gases *************************
346 //-------------------------------------------------
348 //-------------------------------------------------
359 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
369 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
372 //--------------------------------------------------------------
374 //--------------------------------------------------------------
391 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
393 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
408 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
410 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
426 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
428 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
430 //----------------------------------------------------------------
431 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
432 //----------------------------------------------------------------
438 Float_t rho,absl,X0,buf[1];
442 for(nc = 0;nc<fNoComp;nc++)
445 // retrive material constants
447 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
452 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
454 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
455 density += fMixtProp[nc]*rho; // density of the mixture
459 // mixture proportions by weight!
461 for(nc = 0;nc<fNoComp;nc++)
464 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
466 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
467 apure[nnc] : amol[nnc])/am;
471 // Drift gases 1 - nonsensitive, 2 - sensitive
473 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
474 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
483 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
486 //----------------------------------------------------------------------
488 //----------------------------------------------------------------------
510 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
532 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
550 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
568 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
570 // SiO2 - used later for the glass fiber
582 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
591 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
600 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
609 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
627 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
646 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
648 // Epoxy - C14 H20 O3
665 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
673 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
677 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,X0,absl,buf,nbuf);
681 wmat[0]=0.644; // by weight!
684 density=0.5*(1.25+2.265);
686 AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
690 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
692 wmat[0]=0.725; // by weight!
697 AliMixture(39,"G10",amat,zmat,density,2,wmat);
702 //----------------------------------------------------------
703 // tracking media for gases
704 //----------------------------------------------------------
706 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
707 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
708 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
709 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
711 //-----------------------------------------------------------
712 // tracking media for solids
713 //-----------------------------------------------------------
715 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
716 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
717 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
718 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
719 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
720 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
721 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
722 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
723 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
724 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
725 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
726 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
730 void AliTPC::GenerNoise(Int_t tablesize)
733 //Generate table with noise
740 if (fNoiseTable) delete[] fNoiseTable;
741 fNoiseTable = new Float_t[tablesize];
742 fNoiseDepth = tablesize;
743 fCurrentNoise =0; //!index of the noise in the noise table
745 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
746 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
749 Float_t AliTPC::GetNoise()
751 // get noise from table
752 // if ((fCurrentNoise%10)==0)
753 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
754 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
755 return fNoiseTable[fCurrentNoise++];
756 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
760 Bool_t AliTPC::IsSectorActive(Int_t sec)
763 // check if the sector is active
764 if (!fActiveSectors) return kTRUE;
765 else return fActiveSectors[sec];
768 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
770 // activate interesting sectors
771 SetTreeAddress();//just for security
772 if (fActiveSectors) delete [] fActiveSectors;
773 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
774 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
775 for (Int_t i=0;i<n;i++)
776 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
780 void AliTPC::SetActiveSectors(Int_t flag)
783 // activate sectors which were hitted by tracks
785 SetTreeAddress();//just for security
786 if (fHitType==0) return; // if Clones hit - not short volume ID information
787 if (fActiveSectors) delete [] fActiveSectors;
788 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
790 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
793 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
797 Fatal("SetActiveSectors","Can not find TreeH in folder");
800 if (fHitType>1) branch = TreeH()->GetBranch("TPC2");
801 else branch = TreeH()->GetBranch("TPC");
802 Stat_t ntracks = TreeH()->GetEntries();
803 // loop over all hits
804 cout<<"\nAliTPC::SetActiveSectors(): Got "<<ntracks<<" tracks\n";
806 for(Int_t track=0;track<ntracks;track++)
810 if (fTrackHits && fHitType&4) {
811 TBranch * br1 = TreeH()->GetBranch("fVolumes");
812 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
813 br1->GetEvent(track);
814 br2->GetEvent(track);
815 Int_t *volumes = fTrackHits->GetVolumes();
816 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
817 fActiveSectors[volumes[j]]=kTRUE;
821 if (fTrackHitsOld && fHitType&2) {
822 TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
824 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
825 for (UInt_t j=0;j<ar->GetSize();j++){
826 fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
836 void AliTPC::Digits2Clusters(Int_t /*eventnumber*/)
838 //-----------------------------------------------------------------
839 // This is a simple cluster finder.
840 //-----------------------------------------------------------------
841 Error("Digits2Clusters",
842 "Dummy function ! Call AliTPCclusterer::Digits2Clusters(...) instead !");
845 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
846 extern Double_t SigmaZ2(Double_t, Double_t);
847 //_____________________________________________________________________________
848 void AliTPC::Hits2Clusters(Int_t /*eventn*/)
850 //--------------------------------------------------------
851 // TPC simple cluster generator from hits
852 // obtained from the TPC Fast Simulator
853 // The point errors are taken from the parametrization
854 //--------------------------------------------------------
856 //-----------------------------------------------------------------
857 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
858 //-----------------------------------------------------------------
859 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
860 // Jouri.Belikov@cern.ch
861 //----------------------------------------------------------------
863 /////////////////////////////////////////////////////////////////////////////
865 //---------------------------------------------------------------------
866 // ALICE TPC Cluster Parameters
867 //--------------------------------------------------------------------
871 // Cluster width in rphi
872 const Float_t kACrphi=0.18322;
873 const Float_t kBCrphi=0.59551e-3;
874 const Float_t kCCrphi=0.60952e-1;
875 // Cluster width in z
876 const Float_t kACz=0.19081;
877 const Float_t kBCz=0.55938e-3;
878 const Float_t kCCz=0.30428;
882 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
886 //if(fDefaults == 0) SetDefaults();
888 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
890 TParticle *particle; // pointer to a given particle
891 AliTPChit *tpcHit; // pointer to a sigle TPC hit
895 Float_t pl,pt,tanth,rpad,ratio;
898 //---------------------------------------------------------------
899 // Get the access to the tracks
900 //---------------------------------------------------------------
905 Fatal("Hits2Clusters","Can not find TreeH in folder");
910 Stat_t ntracks = tH->GetEntries();
912 //Switch to the output file
914 if (fLoader->TreeR() == 0x0) fLoader->MakeTree("R");
916 cout<<"fTPCParam->GetTitle() = "<<fTPCParam->GetTitle()<<endl;
918 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
920 //fTPCParam->Write(fTPCParam->GetTitle());
922 AliTPCClustersArray carray;
923 carray.Setup(fTPCParam);
924 carray.SetClusterType("AliTPCcluster");
925 carray.MakeTree(fLoader->TreeR());
927 Int_t nclusters=0; //cluster counter
929 //------------------------------------------------------------
930 // Loop over all sectors (72 sectors for 20 deg
931 // segmentation for both lower and upper sectors)
932 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
933 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
935 // First cluster for sector 0 starts at "0"
936 //------------------------------------------------------------
938 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
940 fTPCParam->AdjustCosSin(isec,cph,sph);
942 //------------------------------------------------------------
944 //------------------------------------------------------------
946 for(Int_t track=0;track<ntracks;track++){
950 // Get number of the TPC hits
952 tpcHit = (AliTPChit*)FirstHit(-1);
958 if (tpcHit->fQ == 0.) {
959 tpcHit = (AliTPChit*) NextHit();
960 continue; //information about track (I.Belikov)
962 sector=tpcHit->fSector; // sector number
965 tpcHit = (AliTPChit*) NextHit();
968 ipart=tpcHit->Track();
969 particle=gAlice->GetMCApp()->Particle(ipart);
972 if(pt < 1.e-9) pt=1.e-9;
974 tanth = TMath::Abs(tanth);
975 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
976 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
978 // space-point resolutions
980 sigmaRphi=SigmaY2(rpad,tanth,pt);
981 sigmaZ =SigmaZ2(rpad,tanth );
985 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
986 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
988 // temporary protection
990 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
991 if(sigmaZ < 0.) sigmaZ=0.4e-3;
992 if(clRphi < 0.) clRphi=2.5e-3;
993 if(clZ < 0.) clZ=2.5e-5;
998 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
999 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1001 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
1002 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
1003 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
1004 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
1005 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1006 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
1007 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
1008 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
1009 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
1010 xyz[2]=sigmaRphi; // fSigmaY2
1011 xyz[3]=sigmaZ; // fSigmaZ2
1012 xyz[4]=tpcHit->fQ; // q
1014 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
1015 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
1017 Int_t tracks[3]={tpcHit->Track(), -1, -1};
1018 AliTPCcluster cluster(tracks,xyz);
1020 clrow->InsertCluster(&cluster); nclusters++;
1022 tpcHit = (AliTPChit*)NextHit();
1025 } // end of loop over hits
1027 } // end of loop over tracks
1029 Int_t nrows=fTPCParam->GetNRow(isec);
1030 for (Int_t irow=0; irow<nrows; irow++) {
1031 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
1032 if (!clrow) continue;
1033 carray.StoreRow(isec,irow);
1034 carray.ClearRow(isec,irow);
1037 } // end of loop over sectors
1039 cerr<<"Number of made clusters : "<<nclusters<<" \n";
1040 fLoader->WriteRecPoints("OVERWRITE");
1043 } // end of function
1045 //_________________________________________________________________
1046 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1048 //--------------------------------------------------------
1049 //calculate exact cross point of track and given pad row
1050 //resulting values are expressed in "digit" coordinata
1051 //--------------------------------------------------------
1053 //-----------------------------------------------------------------
1054 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1055 //-----------------------------------------------------------------
1057 if (fClustersArray==0){
1061 TParticle *particle; // pointer to a given particle
1062 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1063 // Int_t sector,nhits;
1065 const Int_t kcmaxhits=30000;
1066 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
1067 AliTPCFastVector & xxx = *xxxx;
1068 Int_t maxhits = kcmaxhits;
1069 //construct array for each padrow
1070 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
1071 fClustersArray->CreateRow(isec,i);
1073 //---------------------------------------------------------------
1074 // Get the access to the tracks
1075 //---------------------------------------------------------------
1077 TTree *tH = TreeH();
1080 Fatal("Hits2Clusters","Can not find TreeH in folder");
1085 Stat_t ntracks = tH->GetEntries();
1086 Int_t npart = gAlice->GetMCApp()->GetNtrack();
1089 if (fHitType>1) branch = tH->GetBranch("TPC2");
1090 else branch = tH->GetBranch("TPC");
1092 //------------------------------------------------------------
1094 //------------------------------------------------------------
1096 for(Int_t track=0;track<ntracks;track++){
1097 Bool_t isInSector=kTRUE;
1099 isInSector = TrackInVolume(isec,track);
1100 if (!isInSector) continue;
1102 branch->GetEntry(track); // get next track
1104 // Get number of the TPC hits and a pointer
1108 Int_t currentIndex=0;
1109 Int_t lastrow=-1; //last writen row
1113 tpcHit = (AliTPChit*)FirstHit(-1);
1116 Int_t sector=tpcHit->fSector; // sector number
1118 tpcHit = (AliTPChit*) NextHit();
1122 ipart=tpcHit->Track();
1123 if (ipart<npart) particle=gAlice->GetMCApp()->Particle(ipart);
1127 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1128 Int_t index[3]={1,isec,0};
1129 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1130 if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
1131 if (lastrow<0) lastrow=currentrow;
1132 if (currentrow==lastrow){
1133 if ( currentIndex>=maxhits){
1135 xxx.ResizeTo(4*maxhits);
1137 xxx(currentIndex*4)=x[0];
1138 xxx(currentIndex*4+1)=x[1];
1139 xxx(currentIndex*4+2)=x[2];
1140 xxx(currentIndex*4+3)=tpcHit->fQ;
1144 if (currentIndex>2){
1156 for (Int_t index=0;index<currentIndex;index++){
1158 x=x2=x3=x4=xxx(index*4);
1166 sumy+=xxx(index*4+1);
1167 sumxy+=xxx(index*4+1)*x;
1168 sumx2y+=xxx(index*4+1)*x2;
1169 sumz+=xxx(index*4+2);
1170 sumxz+=xxx(index*4+2)*x;
1171 sumx2z+=xxx(index*4+2)*x2;
1172 sumq+=xxx(index*4+3);
1174 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1175 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1176 sumx2*(sumx*sumx3-sumx2*sumx2);
1178 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1179 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1180 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1181 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1183 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1184 sumx2*(sumx*sumx2y-sumx2*sumxy);
1185 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1186 sumx2*(sumx*sumx2z-sumx2*sumxz);
1188 if (TMath::Abs(det)<0.00001){
1189 tpcHit = (AliTPChit*)NextHit();
1193 Float_t y=detay/det+centralPad;
1194 Float_t z=detaz/det;
1195 Float_t by=detby/det; //y angle
1196 Float_t bz=detbz/det; //z angle
1197 sumy/=Float_t(currentIndex);
1198 sumz/=Float_t(currentIndex);
1200 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1202 AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
1203 // AliComplexCluster cl;
1209 cl->fTracks[0]=ipart;
1213 } //end of calculating cluster for given row
1216 tpcHit = (AliTPChit*)NextHit();
1217 } // end of loop over hits
1218 } // end of loop over tracks
1219 //write padrows to tree
1220 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1221 fClustersArray->StoreRow(isec,ii);
1222 fClustersArray->ClearRow(isec,ii);
1231 void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
1233 //create digits from summable digits
1234 GenerNoise(500000); //create teble with noise
1236 //conect tree with sSDigits
1237 TTree *t = fLoader->TreeS();
1241 fLoader->LoadSDigits("READ");
1242 t = fLoader->TreeS();
1245 Error("SDigits2Digits2","Can not get input TreeS");
1250 if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
1252 AliSimDigits digarr, *dummy=&digarr;
1253 TBranch* sdb = t->GetBranch("Segment");
1256 Error("SDigits2Digits2","Can not find branch with segments in TreeS.");
1260 sdb->SetAddress(&dummy);
1262 Stat_t nentries = t->GetEntries();
1264 // set zero suppression
1266 fTPCParam->SetZeroSup(2);
1268 // get zero suppression
1270 Int_t zerosup = fTPCParam->GetZeroSup();
1272 //make tree with digits
1274 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1275 arr->SetClass("AliSimDigits");
1276 arr->Setup(fTPCParam);
1277 arr->MakeTree(fLoader->TreeD());
1279 AliTPCParam * par = fTPCParam;
1281 //Loop over segments of the TPC
1283 for (Int_t n=0; n<nentries; n++) {
1286 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1287 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
1290 if (!IsSectorActive(sec))
1292 // cout<<n<<" NOT Active \n";
1297 // cout<<n<<" Active \n";
1299 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1300 Int_t nrows = digrow->GetNRows();
1301 Int_t ncols = digrow->GetNCols();
1303 digrow->ExpandBuffer();
1304 digarr.ExpandBuffer();
1305 digrow->ExpandTrackBuffer();
1306 digarr.ExpandTrackBuffer();
1309 Short_t * pamp0 = digarr.GetDigits();
1310 Int_t * ptracks0 = digarr.GetTracks();
1311 Short_t * pamp1 = digrow->GetDigits();
1312 Int_t * ptracks1 = digrow->GetTracks();
1313 Int_t nelems =nrows*ncols;
1314 Int_t saturation = fTPCParam->GetADCSat();
1315 //use internal structure of the AliDigits - for speed reason
1316 //if you cahnge implementation
1317 //of the Alidigits - it must be rewriten -
1318 for (Int_t i= 0; i<nelems; i++){
1319 // Float_t q = *pamp0;
1320 //q/=16.; //conversion faktor
1321 //Float_t noise= GetNoise();
1323 //q= TMath::Nint(q);
1324 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1326 if (q>saturation) q=saturation;
1328 //if (ptracks0[0]==0)
1331 ptracks1[0]=ptracks0[0];
1332 ptracks1[nelems]=ptracks0[nelems];
1333 ptracks1[2*nelems]=ptracks0[2*nelems];
1341 arr->StoreRow(sec,row);
1342 arr->ClearRow(sec,row);
1343 // cerr<<sec<<"\t"<<row<<"\n";
1348 fLoader->WriteDigits("OVERWRITE");
1352 //__________________________________________________________________
1353 void AliTPC::SetDefaults(){
1356 cerr<<"Setting default parameters...\n";
1358 // Set response functions
1361 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1363 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1365 printf("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...\n");
1367 param = new AliTPCParamSR();
1370 param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1373 printf("No TPC parameters found\n");
1378 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1379 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1380 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1381 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1382 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1383 rf->SetOffset(3*param->GetZSigma());
1386 TDirectory *savedir=gDirectory;
1387 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1389 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1392 prfinner->Read("prf_07504_Gati_056068_d02");
1393 prfouter1->Read("prf_10006_Gati_047051_d03");
1394 prfouter2->Read("prf_15006_Gati_047051_d03");
1398 param->SetInnerPRF(prfinner);
1399 param->SetOuter1PRF(prfouter1);
1400 param->SetOuter2PRF(prfouter2);
1401 param->SetTimeRF(rf);
1411 //__________________________________________________________________
1412 void AliTPC::Hits2Digits(Int_t eventnumber)
1414 //----------------------------------------------------
1415 // Loop over all sectors for a single event
1416 //----------------------------------------------------
1417 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1418 rl->GetEvent(eventnumber);
1419 if (fLoader->TreeH() == 0x0)
1421 if(fLoader->LoadHits())
1423 Error("Hits2Digits","Can not load hits.");
1428 if (fLoader->TreeD() == 0x0 )
1430 fLoader->MakeTree("D");
1431 if (fLoader->TreeD() == 0x0 )
1433 Error("Hits2Digits","Can not get TreeD");
1438 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1439 GenerNoise(500000); //create teble with noise
1441 //setup TPCDigitsArray
1443 if(GetDigitsArray()) delete GetDigitsArray();
1445 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1446 arr->SetClass("AliSimDigits");
1447 arr->Setup(fTPCParam);
1449 arr->MakeTree(fLoader->TreeD());
1450 SetDigitsArray(arr);
1452 fDigitsSwitch=0; // standard digits
1454 cerr<<"Digitizing TPC -- normal digits...\n";
1456 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1457 if (IsSectorActive(isec))
1459 if (fDebug) Info("Hits2Digits","Sector %d is active.",isec);
1460 Hits2DigitsSector(isec);
1464 if (fDebug) Info("Hits2Digits","Sector %d is NOT active.",isec);
1467 fLoader->WriteDigits("OVERWRITE");
1469 //this line prevents the crash in the similar one
1470 //on the beginning of this method
1471 //destructor attempts to reset the tree, which is deleted by the loader
1472 //need to be redesign
1473 if(GetDigitsArray()) delete GetDigitsArray();
1474 SetDigitsArray(0x0);
1478 //__________________________________________________________________
1479 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1482 //-----------------------------------------------------------
1483 // summable digits - 16 bit "ADC", no noise, no saturation
1484 //-----------------------------------------------------------
1486 //----------------------------------------------------
1487 // Loop over all sectors for a single event
1488 //----------------------------------------------------
1489 // AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1491 AliRunLoader* rl = fLoader->GetRunLoader();
1493 rl->GetEvent(eventnumber);
1494 if (fLoader->TreeH() == 0x0)
1496 if(fLoader->LoadHits())
1498 Error("Hits2Digits","Can not load hits.");
1505 if (fLoader->TreeS() == 0x0 )
1507 fLoader->MakeTree("S");
1510 if(fDefaults == 0) SetDefaults();
1512 GenerNoise(500000); //create table with noise
1513 //setup TPCDigitsArray
1515 if(GetDigitsArray()) delete GetDigitsArray();
1518 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1519 arr->SetClass("AliSimDigits");
1520 arr->Setup(fTPCParam);
1521 arr->MakeTree(fLoader->TreeS());
1523 SetDigitsArray(arr);
1525 cerr<<"Digitizing TPC -- summable digits...\n";
1527 fDigitsSwitch=1; // summable digits
1529 // set zero suppression to "0"
1531 fTPCParam->SetZeroSup(0);
1533 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1534 if (IsSectorActive(isec))
1536 // cout<<"Sector "<<isec<<" is active\n";
1537 Hits2DigitsSector(isec);
1540 fLoader->WriteSDigits("OVERWRITE");
1542 //this line prevents the crash in the similar one
1543 //on the beginning of this method
1544 //destructor attempts to reset the tree, which is deleted by the loader
1545 //need to be redesign
1546 if(GetDigitsArray()) delete GetDigitsArray();
1547 SetDigitsArray(0x0);
1549 //__________________________________________________________________
1551 void AliTPC::Hits2SDigits()
1554 //-----------------------------------------------------------
1555 // summable digits - 16 bit "ADC", no noise, no saturation
1556 //-----------------------------------------------------------
1558 //----------------------------------------------------
1559 // Loop over all sectors for a single event
1560 //----------------------------------------------------
1561 //MI change - for pp run
1562 // Int_t eventnumber = gAlice->GetEvNumber();
1563 // AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1564 // rl->GetEvent(eventnumber);
1565 // 12/05/2003 This method is obsolete and not used. It should be redesingned
1568 if (fLoader->TreeH() == 0x0)
1570 if(fLoader->LoadHits())
1572 Error("Hits2Digits","Can not load hits.");
1577 if(fDefaults == 0) SetDefaults();
1578 GenerNoise(500000); //create table with noise
1580 //setup TPCDigitsArray
1582 if(GetDigitsArray()) delete GetDigitsArray();
1584 if (fLoader->TreeS() == 0x0 )
1586 fLoader->MakeTree("S");
1589 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1590 arr->SetClass("AliSimDigits");
1591 arr->Setup(fTPCParam);
1592 arr->MakeTree(fLoader->TreeS());
1593 SetDigitsArray(arr);
1595 cerr<<"Digitizing TPC -- summable digits...\n";
1597 // fDigitsSwitch=1; // summable digits -for the moment direct
1599 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1603 cout<<"Why method TPC::Hits2SDigits writes digits and not sdigits? skowron\n";
1604 fLoader->WriteDigits("OVERWRITE");
1606 //_____________________________________________________________________________
1608 void AliTPC::Hits2DigitsSector(Int_t isec)
1610 //-------------------------------------------------------------------
1611 // TPC conversion from hits to digits.
1612 //-------------------------------------------------------------------
1614 //-----------------------------------------------------------------
1615 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1616 //-----------------------------------------------------------------
1618 //-------------------------------------------------------
1619 // Get the access to the track hits
1620 //-------------------------------------------------------
1622 // check if the parameters are set - important if one calls this method
1623 // directly, not from the Hits2Digits
1625 if(fDefaults == 0) SetDefaults();
1627 TTree *tH = TreeH(); // pointer to the hits tree
1630 Fatal("Hits2DigitsSector","Can not find TreeH in folder");
1634 Stat_t ntracks = tH->GetEntries();
1638 //-------------------------------------------
1639 // Only if there are any tracks...
1640 //-------------------------------------------
1644 //printf("*** Processing sector number %d ***\n",isec);
1646 Int_t nrows =fTPCParam->GetNRow(isec);
1648 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1650 MakeSector(isec,nrows,tH,ntracks,row);
1652 //--------------------------------------------------------
1653 // Digitize this sector, row by row
1654 // row[i] is the pointer to the TObjArray of AliTPCFastVectors,
1655 // each one containing electrons accepted on this
1656 // row, assigned into tracks
1657 //--------------------------------------------------------
1661 if (fDigitsArray->GetTree()==0)
1663 Fatal("Hits2DigitsSector","Tree not set in fDigitsArray");
1666 for (i=0;i<nrows;i++){
1668 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1670 DigitizeRow(i,isec,row);
1672 fDigitsArray->StoreRow(isec,i);
1674 Int_t ndig = dig->GetDigitSize();
1677 printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1679 fDigitsArray->ClearRow(isec,i);
1682 } // end of the sector digitization
1684 for(i=0;i<nrows+2;i++){
1689 delete [] row; // delete the array of pointers to TObjArray-s
1693 } // end of Hits2DigitsSector
1696 //_____________________________________________________________________________
1697 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1699 //-----------------------------------------------------------
1700 // Single row digitization, coupling from the neighbouring
1701 // rows taken into account
1702 //-----------------------------------------------------------
1704 //-----------------------------------------------------------------
1705 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1706 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1707 //-----------------------------------------------------------------
1710 Float_t zerosup = fTPCParam->GetZeroSup();
1711 // Int_t nrows =fTPCParam->GetNRow(isec);
1712 fCurrentIndex[1]= isec;
1715 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1716 Int_t nofTbins = fTPCParam->GetMaxTBin();
1717 Int_t indexRange[4];
1719 // Integrated signal for this row
1720 // and a single track signal
1723 AliTPCFastMatrix *m1 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // integrated
1724 AliTPCFastMatrix *m2 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // single
1726 AliTPCFastMatrix &total = *m1;
1728 // Array of pointers to the label-signal list
1730 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1731 Float_t **pList = new Float_t* [nofDigits];
1735 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1739 //Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1740 //Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1743 for (Int_t row= row1;row<=row2;row++){
1744 Int_t nTracks= rows[row]->GetEntries();
1745 for (i1=0;i1<nTracks;i1++){
1746 fCurrentIndex[2]= row;
1747 fCurrentIndex[3]=irow+1;
1749 m2->Zero(); // clear single track signal matrix
1750 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1751 GetList(trackLabel,nofPads,m2,indexRange,pList);
1753 else GetSignal(rows[row],i1,0,m1,indexRange);
1759 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1761 Float_t fzerosup = zerosup+0.5;
1762 for(Int_t it=0;it<nofTbins;it++){
1763 Float_t *pq = &(total.UncheckedAt(0,it));
1764 for(Int_t ip=0;ip<nofPads;ip++){
1768 if(fDigitsSwitch == 0){
1770 if(q <=fzerosup) continue; // do not fill zeros
1772 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1777 if(q <= 0.) continue; // do not fill zeros
1778 if(q>2000.) q=2000.;
1784 // "real" signal or electronic noise (list = -1)?
1787 for(Int_t j1=0;j1<3;j1++){
1788 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1793 <A NAME="AliDigits"></A>
1794 using of AliDigits object
1797 dig->SetDigitFast((Short_t)q,it,ip);
1798 if (fDigitsArray->IsSimulated())
1800 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1801 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1802 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1806 } // end of loop over time buckets
1807 } // end of lop over pads
1810 // This row has been digitized, delete nonused stuff
1813 for(lp=0;lp<nofDigits;lp++){
1814 if(pList[lp]) delete [] pList[lp];
1823 } // end of DigitizeRow
1825 //_____________________________________________________________________________
1827 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1828 AliTPCFastMatrix *m1, AliTPCFastMatrix *m2,Int_t *indexRange)
1831 //---------------------------------------------------------------
1832 // Calculates 2-D signal (pad,time) for a single track,
1833 // returns a pointer to the signal matrix and the track label
1834 // No digitization is performed at this level!!!
1835 //---------------------------------------------------------------
1837 //-----------------------------------------------------------------
1838 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1839 // Modified: Marian Ivanov
1840 //-----------------------------------------------------------------
1842 AliTPCFastVector *tv;
1844 tv = (AliTPCFastVector*)p1->At(ntr); // pointer to a track
1845 AliTPCFastVector &v = *tv;
1847 Float_t label = v(0);
1848 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1)-1)/2;
1850 Int_t nElectrons = (tv->GetNrows()-1)/4;
1851 indexRange[0]=9999; // min pad
1852 indexRange[1]=-1; // max pad
1853 indexRange[2]=9999; //min time
1854 indexRange[3]=-1; // max time
1856 AliTPCFastMatrix &signal = *m1;
1857 AliTPCFastMatrix &total = *m2;
1859 // Loop over all electrons
1861 for(Int_t nel=0; nel<nElectrons; nel++){
1863 Float_t aval = v(idx+4);
1864 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1865 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1866 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
1868 Int_t *index = fTPCParam->GetResBin(0);
1869 Float_t *weight = & (fTPCParam->GetResWeight(0));
1871 if (n>0) for (Int_t i =0; i<n; i++){
1872 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1875 Int_t time=index[2];
1876 Float_t qweight = *(weight)*eltoadcfac;
1878 if (m1!=0) signal.UncheckedAt(pad,time)+=qweight;
1879 total.UncheckedAt(pad,time)+=qweight;
1880 if (indexRange[0]>pad) indexRange[0]=pad;
1881 if (indexRange[1]<pad) indexRange[1]=pad;
1882 if (indexRange[2]>time) indexRange[2]=time;
1883 if (indexRange[3]<time) indexRange[3]=time;
1890 } // end of loop over electrons
1892 return label; // returns track label when finished
1895 //_____________________________________________________________________________
1896 void AliTPC::GetList(Float_t label,Int_t np,AliTPCFastMatrix *m,
1897 Int_t *indexRange, Float_t **pList)
1899 //----------------------------------------------------------------------
1900 // Updates the list of tracks contributing to digits for a given row
1901 //----------------------------------------------------------------------
1903 //-----------------------------------------------------------------
1904 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1905 //-----------------------------------------------------------------
1907 AliTPCFastMatrix &signal = *m;
1909 // lop over nonzero digits
1911 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1912 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1915 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1917 if(signal(ip,it)<0.5) continue;
1920 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1922 if(!pList[globalIndex]){
1925 // Create new list (6 elements - 3 signals and 3 labels),
1928 pList[globalIndex] = new Float_t [6];
1932 *pList[globalIndex] = -1.;
1933 *(pList[globalIndex]+1) = -1.;
1934 *(pList[globalIndex]+2) = -1.;
1935 *(pList[globalIndex]+3) = -1.;
1936 *(pList[globalIndex]+4) = -1.;
1937 *(pList[globalIndex]+5) = -1.;
1940 *pList[globalIndex] = label;
1941 *(pList[globalIndex]+3) = signal(ip,it);
1945 // check the signal magnitude
1947 Float_t highest = *(pList[globalIndex]+3);
1948 Float_t middle = *(pList[globalIndex]+4);
1949 Float_t lowest = *(pList[globalIndex]+5);
1952 // compare the new signal with already existing list
1955 if(signal(ip,it)<lowest) continue; // neglect this track
1959 if (signal(ip,it)>highest){
1960 *(pList[globalIndex]+5) = middle;
1961 *(pList[globalIndex]+4) = highest;
1962 *(pList[globalIndex]+3) = signal(ip,it);
1964 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1965 *(pList[globalIndex]+1) = *pList[globalIndex];
1966 *pList[globalIndex] = label;
1968 else if (signal(ip,it)>middle){
1969 *(pList[globalIndex]+5) = middle;
1970 *(pList[globalIndex]+4) = signal(ip,it);
1972 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
1973 *(pList[globalIndex]+1) = label;
1976 *(pList[globalIndex]+5) = signal(ip,it);
1977 *(pList[globalIndex]+2) = label;
1981 } // end of loop over pads
1982 } // end of loop over time bins
1987 //___________________________________________________________________
1988 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
1989 Stat_t ntracks,TObjArray **row)
1992 //-----------------------------------------------------------------
1993 // Prepares the sector digitization, creates the vectors of
1994 // tracks for each row of this sector. The track vector
1995 // contains the track label and the position of electrons.
1996 //-----------------------------------------------------------------
1998 //-----------------------------------------------------------------
1999 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2000 //-----------------------------------------------------------------
2002 Float_t gasgain = fTPCParam->GetGasGain();
2006 AliTPChit *tpcHit; // pointer to a sigle TPC hit
2009 if (fHitType>1) branch = TH->GetBranch("TPC2");
2010 else branch = TH->GetBranch("TPC");
2013 //----------------------------------------------
2014 // Create TObjArray-s, one for each row,
2015 // each TObjArray will store the AliTPCFastVectors
2016 // of electrons, one AliTPCFastVectors per each track.
2017 //----------------------------------------------
2019 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
2020 AliTPCFastVector **tracks = new AliTPCFastVector* [nrows+2]; //pointers to the track vectors
2022 for(i=0; i<nrows+2; i++){
2023 row[i] = new TObjArray;
2030 //--------------------------------------------------------------------
2031 // Loop over tracks, the "track" contains the full history
2032 //--------------------------------------------------------------------
2034 Int_t previousTrack,currentTrack;
2035 previousTrack = -1; // nothing to store so far!
2037 for(Int_t track=0;track<ntracks;track++){
2038 Bool_t isInSector=kTRUE;
2040 isInSector = TrackInVolume(isec,track);
2041 if (!isInSector) continue;
2043 branch->GetEntry(track); // get next track
2047 tpcHit = (AliTPChit*)FirstHit(-1);
2049 //--------------------------------------------------------------
2051 //--------------------------------------------------------------
2056 Int_t sector=tpcHit->fSector; // sector number
2058 tpcHit = (AliTPChit*) NextHit();
2062 currentTrack = tpcHit->Track(); // track number
2065 if(currentTrack != previousTrack){
2067 // store already filled fTrack
2069 for(i=0;i<nrows+2;i++){
2070 if(previousTrack != -1){
2071 if(nofElectrons[i]>0){
2072 AliTPCFastVector &v = *tracks[i];
2073 v(0) = previousTrack;
2074 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2075 row[i]->Add(tracks[i]);
2078 delete tracks[i]; // delete empty AliTPCFastVector
2084 tracks[i] = new AliTPCFastVector(481); // AliTPCFastVectors for the next fTrack
2086 } // end of loop over rows
2088 previousTrack=currentTrack; // update track label
2091 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2093 //---------------------------------------------------
2094 // Calculate the electron attachment probability
2095 //---------------------------------------------------
2098 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
2099 /fTPCParam->GetDriftV();
2101 Float_t attProb = fTPCParam->GetAttCoef()*
2102 fTPCParam->GetOxyCont()*time; // fraction!
2104 //-----------------------------------------------
2105 // Loop over electrons
2106 //-----------------------------------------------
2109 for(Int_t nel=0;nel<qI;nel++){
2110 // skip if electron lost due to the attachment
2111 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2116 // protection for the nonphysical avalanche size (10**6 maximum)
2118 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2119 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2122 TransportElectron(xyz,index);
2124 fTPCParam->GetPadRow(xyz,index);
2125 // row 0 - cross talk from the innermost row
2126 // row fNRow+1 cross talk from the outermost row
2127 rowNumber = index[2]+1;
2128 //transform position to local digit coordinates
2129 //relative to nearest pad row
2130 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2132 if (isec <fTPCParam->GetNInnerSector()) {
2133 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2134 y1 = fTPCParam->GetYInner(rowNumber);
2137 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2138 y1 = fTPCParam->GetYOuter(rowNumber);
2140 // gain inefficiency at the wires edges - linear
2143 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));
2145 nofElectrons[rowNumber]++;
2146 //----------------------------------
2147 // Expand vector if necessary
2148 //----------------------------------
2149 if(nofElectrons[rowNumber]>120){
2150 Int_t range = tracks[rowNumber]->GetNrows();
2151 if((nofElectrons[rowNumber])>(range-1)/4){
2153 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
2157 AliTPCFastVector &v = *tracks[rowNumber];
2158 Int_t idx = 4*nofElectrons[rowNumber]-3;
2159 Real_t * position = &(((AliTPCFastVector&)v).UncheckedAt(idx)); //make code faster
2160 memcpy(position,xyz,4*sizeof(Float_t));
2162 } // end of loop over electrons
2164 tpcHit = (AliTPChit*)NextHit();
2166 } // end of loop over hits
2167 } // end of loop over tracks
2170 // store remaining track (the last one) if not empty
2173 for(i=0;i<nrows+2;i++){
2174 if(nofElectrons[i]>0){
2175 AliTPCFastVector &v = *tracks[i];
2176 v(0) = previousTrack;
2177 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2178 row[i]->Add(tracks[i]);
2187 delete [] nofElectrons;
2190 } // end of MakeSector
2193 //_____________________________________________________________________________
2197 // Initialise TPC detector after definition of geometry
2202 printf("\n%s: ",ClassName());
2203 for(i=0;i<35;i++) printf("*");
2204 printf(" TPC_INIT ");
2205 for(i=0;i<35;i++) printf("*");
2206 printf("\n%s: ",ClassName());
2208 for(i=0;i<80;i++) printf("*");
2213 //_____________________________________________________________________________
2214 void AliTPC::MakeBranch(Option_t* option)
2217 // Create Tree branches for the TPC.
2219 if(GetDebug()) Info("MakeBranch","");
2220 Int_t buffersize = 4000;
2221 char branchname[10];
2222 sprintf(branchname,"%s",GetName());
2224 const char *h = strstr(option,"H");
2226 if ( h && (fHitType<=1) && (fHits == 0x0)) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2228 AliDetector::MakeBranch(option);
2230 const char *d = strstr(option,"D");
2232 if (fDigits && fLoader->TreeD() && d)
2234 MakeBranchInTree(gAlice->TreeD(), branchname, &fDigits, buffersize, 0);
2237 if (fHitType>1) MakeBranch2(option,0); // MI change 14.09.2000
2240 //_____________________________________________________________________________
2241 void AliTPC::ResetDigits()
2244 // Reset number of digits and the digits array for this detector
2247 if (fDigits) fDigits->Clear();
2250 //_____________________________________________________________________________
2251 void AliTPC::SetSecAL(Int_t sec)
2253 //---------------------------------------------------
2254 // Activate/deactivate selection for lower sectors
2255 //---------------------------------------------------
2257 //-----------------------------------------------------------------
2258 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2259 //-----------------------------------------------------------------
2263 //_____________________________________________________________________________
2264 void AliTPC::SetSecAU(Int_t sec)
2266 //----------------------------------------------------
2267 // Activate/deactivate selection for upper sectors
2268 //---------------------------------------------------
2270 //-----------------------------------------------------------------
2271 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2272 //-----------------------------------------------------------------
2276 //_____________________________________________________________________________
2277 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2279 //----------------------------------------
2280 // Select active lower sectors
2281 //----------------------------------------
2283 //-----------------------------------------------------------------
2284 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2285 //-----------------------------------------------------------------
2295 //_____________________________________________________________________________
2296 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2297 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2298 Int_t s11 , Int_t s12)
2300 //--------------------------------
2301 // Select active upper sectors
2302 //--------------------------------
2304 //-----------------------------------------------------------------
2305 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2306 //-----------------------------------------------------------------
2322 //_____________________________________________________________________________
2323 void AliTPC::SetSens(Int_t sens)
2326 //-------------------------------------------------------------
2327 // Activates/deactivates the sensitive strips at the center of
2328 // the pad row -- this is for the space-point resolution calculations
2329 //-------------------------------------------------------------
2331 //-----------------------------------------------------------------
2332 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2333 //-----------------------------------------------------------------
2339 void AliTPC::SetSide(Float_t side=0.)
2341 // choice of the TPC side
2346 //____________________________________________________________________________
2347 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2348 Float_t p2,Float_t p3)
2351 // gax mixture definition
2365 //_____________________________________________________________________________
2367 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2370 // electron transport taking into account:
2372 // 2.ExB at the wires
2373 // 3. nonisochronity
2375 // xyz and index must be already transformed to system 1
2378 fTPCParam->Transform1to2(xyz,index);
2381 Float_t driftl=xyz[2];
2382 if(driftl<0.01) driftl=0.01;
2383 driftl=TMath::Sqrt(driftl);
2384 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2385 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2386 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2387 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2388 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2392 if (fTPCParam->GetMWPCReadout()==kTRUE){
2393 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2394 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2396 //add nonisochronity (not implemented yet)
2399 ClassImp(AliTPCdigit)
2401 //_____________________________________________________________________________
2402 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2406 // Creates a TPC digit object
2408 fSector = digits[0];
2409 fPadRow = digits[1];
2412 fSignal = digits[4];
2418 //_____________________________________________________________________________
2419 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2423 // Creates a TPC hit object
2433 //________________________________________________________________________
2434 // Additional code because of the AliTPCTrackHitsV2
2436 void AliTPC::MakeBranch2(Option_t *option,const char */*file*/)
2439 // Create a new branch in the current Root Tree
2440 // The branch of fHits is automatically split
2441 // MI change 14.09.2000
2442 if(GetDebug()) Info("MakeBranch2","");
2443 if (fHitType<2) return;
2444 char branchname[10];
2445 sprintf(branchname,"%s2",GetName());
2447 // Get the pointer to the header
2448 const char *cH = strstr(option,"H");
2450 if (fTrackHits && TreeH() && cH && fHitType&4)
2452 if(GetDebug()) Info("MakeBranch2","Making branch for Type 4 Hits");
2453 TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2456 if (fTrackHitsOld && TreeH() && cH && fHitType&2)
2458 if(GetDebug()) Info("MakeBranch2","Making branch for Type 2 Hits");
2459 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2460 TreeH(),fBufferSize,99);
2461 TreeH()->GetListOfBranches()->Add(branch);
2465 void AliTPC::SetTreeAddress()
2467 //Sets tree address for hits
2470 if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2471 AliDetector::SetTreeAddress();
2473 if (fHitType>1) SetTreeAddress2();
2476 void AliTPC::SetTreeAddress2()
2479 // Set branch address for the TrackHits Tree
2481 if(GetDebug()) Info("SetTreeAddress2","");
2484 char branchname[20];
2485 sprintf(branchname,"%s2",GetName());
2487 // Branch address for hit tree
2488 TTree *treeH = TreeH();
2489 if ((treeH)&&(fHitType&4)) {
2490 branch = treeH->GetBranch(branchname);
2493 branch->SetAddress(&fTrackHits);
2494 if (GetDebug()) Info("SetTreeAddress2","fHitType&4 Setting");
2497 if (GetDebug()) Info("SetTreeAddress2","fHitType&4 Failed (can not find branch)");
2500 if ((treeH)&&(fHitType&2)) {
2501 branch = treeH->GetBranch(branchname);
2504 branch->SetAddress(&fTrackHitsOld);
2505 if (GetDebug()) Info("SetTreeAddress2","fHitType&2 Setting");
2507 else if (GetDebug())
2508 Info("SetTreeAddress2","fHitType&2 Failed (can not find branch)");
2510 //set address to TREETR
2512 TTree *treeTR = TreeTR();
2513 if (treeTR && fTrackReferences) {
2514 branch = treeTR->GetBranch(GetName());
2515 if (branch) branch->SetAddress(&fTrackReferences);
2520 void AliTPC::FinishPrimary()
2522 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2523 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2527 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2530 // add hit to the list
2533 int primary = gAlice->GetMCApp()->GetPrimary(track);
2534 gAlice->GetMCApp()->Particle(primary)->SetBit(kKeepBit);
2538 gAlice->GetMCApp()->FlagTrack(track);
2540 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2541 //if (hit->fTrack!=rtrack)
2542 // cout<<"bad track number\n";
2543 if (fTrackHits && fHitType&4)
2544 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2545 hits[1],hits[2],(Int_t)hits[3]);
2546 if (fTrackHitsOld &&fHitType&2 )
2547 fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2548 hits[1],hits[2],(Int_t)hits[3]);
2552 void AliTPC::ResetHits()
2554 if (fHitType&1) AliDetector::ResetHits();
2555 if (fHitType>1) ResetHits2();
2558 void AliTPC::ResetHits2()
2562 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2563 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2567 AliHit* AliTPC::FirstHit(Int_t track)
2569 if (fHitType>1) return FirstHit2(track);
2570 return AliDetector::FirstHit(track);
2572 AliHit* AliTPC::NextHit()
2574 if (fHitType>1) return NextHit2();
2576 return AliDetector::NextHit();
2579 AliHit* AliTPC::FirstHit2(Int_t track)
2582 // Initialise the hit iterator
2583 // Return the address of the first hit for track
2584 // If track>=0 the track is read from disk
2585 // while if track<0 the first hit of the current
2586 // track is returned
2589 gAlice->ResetHits();
2590 TreeH()->GetEvent(track);
2593 if (fTrackHits && fHitType&4) {
2594 fTrackHits->First();
2595 return fTrackHits->GetHit();
2597 if (fTrackHitsOld && fHitType&2) {
2598 fTrackHitsOld->First();
2599 return fTrackHitsOld->GetHit();
2605 AliHit* AliTPC::NextHit2()
2608 //Return the next hit for the current track
2611 if (fTrackHitsOld && fHitType&2) {
2612 fTrackHitsOld->Next();
2613 return fTrackHitsOld->GetHit();
2617 return fTrackHits->GetHit();
2623 void AliTPC::LoadPoints(Int_t)
2627 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2630 if(fHitType==1) AliDetector::LoadPoints(a);
2631 else LoadPoints2(a);
2638 void AliTPC::RemapTrackHitIDs(Int_t *map)
2640 if (!fTrackHits) return;
2642 if (fTrackHitsOld && fHitType&2){
2643 AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2644 for (UInt_t i=0;i<arr->GetSize();i++){
2645 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2646 info->fTrackID = map[info->fTrackID];
2649 if (fTrackHitsOld && fHitType&4){
2650 TClonesArray * arr = fTrackHits->GetArray();;
2651 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2652 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2653 info->fTrackID = map[info->fTrackID];
2658 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2660 //return bool information - is track in given volume
2661 //load only part of the track information
2662 //return true if current track is in volume
2665 if (fTrackHitsOld && fHitType&2) {
2666 TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
2667 br->GetEvent(track);
2668 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2669 for (UInt_t j=0;j<ar->GetSize();j++){
2670 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2674 if (fTrackHits && fHitType&4) {
2675 TBranch * br1 = TreeH()->GetBranch("fVolumes");
2676 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
2677 br2->GetEvent(track);
2678 br1->GetEvent(track);
2679 Int_t *volumes = fTrackHits->GetVolumes();
2680 Int_t nvolumes = fTrackHits->GetNVolumes();
2681 if (!volumes && nvolumes>0) {
2682 printf("Problematic track\t%d\t%d",track,nvolumes);
2685 for (Int_t j=0;j<nvolumes; j++)
2686 if (volumes[j]==id) return kTRUE;
2690 TBranch * br = TreeH()->GetBranch("fSector");
2691 br->GetEvent(track);
2692 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2693 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2700 //_____________________________________________________________________________
2701 void AliTPC::LoadPoints2(Int_t)
2704 // Store x, y, z of all hits in memory
2706 if (fTrackHits == 0 && fTrackHitsOld==0) return;
2709 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2710 if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2712 if (nhits == 0) return;
2713 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2714 if (fPoints == 0) fPoints = new TObjArray(tracks);
2717 Int_t *ntrk=new Int_t[tracks];
2718 Int_t *limi=new Int_t[tracks];
2719 Float_t **coor=new Float_t*[tracks];
2720 for(Int_t i=0;i<tracks;i++) {
2726 AliPoints *points = 0;
2729 Int_t chunk=nhits/4+1;
2731 // Loop over all the hits and store their position
2733 ahit = FirstHit2(-1);
2735 trk=ahit->GetTrack();
2736 if(ntrk[trk]==limi[trk]) {
2738 // Initialise a new track
2739 fp=new Float_t[3*(limi[trk]+chunk)];
2741 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2742 delete [] coor[trk];
2749 fp[3*ntrk[trk] ] = ahit->X();
2750 fp[3*ntrk[trk]+1] = ahit->Y();
2751 fp[3*ntrk[trk]+2] = ahit->Z();
2759 for(trk=0; trk<tracks; ++trk) {
2761 points = new AliPoints();
2762 points->SetMarkerColor(GetMarkerColor());
2763 points->SetMarkerSize(GetMarkerSize());
2764 points->SetDetector(this);
2765 points->SetParticle(trk);
2766 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2767 fPoints->AddAt(points,trk);
2768 delete [] coor[trk];
2778 //_____________________________________________________________________________
2779 void AliTPC::LoadPoints3(Int_t)
2782 // Store x, y, z of all hits in memory
2783 // - only intersection point with pad row
2784 if (fTrackHits == 0) return;
2786 Int_t nhits = fTrackHits->GetEntriesFast();
2787 if (nhits == 0) return;
2788 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2789 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2790 fPoints->Expand(2*tracks);
2793 Int_t *ntrk=new Int_t[tracks];
2794 Int_t *limi=new Int_t[tracks];
2795 Float_t **coor=new Float_t*[tracks];
2796 for(Int_t i=0;i<tracks;i++) {
2802 AliPoints *points = 0;
2805 Int_t chunk=nhits/4+1;
2807 // Loop over all the hits and store their position
2809 ahit = FirstHit2(-1);
2810 //for (Int_t hit=0;hit<nhits;hit++) {
2814 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2815 trk=ahit->GetTrack();
2816 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2817 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2818 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2819 if (currentrow!=lastrow){
2820 lastrow = currentrow;
2821 //later calculate intersection point
2822 if(ntrk[trk]==limi[trk]) {
2824 // Initialise a new track
2825 fp=new Float_t[3*(limi[trk]+chunk)];
2827 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2828 delete [] coor[trk];
2835 fp[3*ntrk[trk] ] = ahit->X();
2836 fp[3*ntrk[trk]+1] = ahit->Y();
2837 fp[3*ntrk[trk]+2] = ahit->Z();
2844 for(trk=0; trk<tracks; ++trk) {
2846 points = new AliPoints();
2847 points->SetMarkerColor(GetMarkerColor()+1);
2848 points->SetMarkerStyle(5);
2849 points->SetMarkerSize(0.2);
2850 points->SetDetector(this);
2851 points->SetParticle(trk);
2852 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2853 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2854 fPoints->AddAt(points,tracks+trk);
2855 delete [] coor[trk];
2866 void AliTPC::FindTrackHitsIntersection(TClonesArray * /*arr*/)
2870 //fill clones array with intersection of current point with the
2875 const Int_t kcmaxhits=30000;
2876 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
2877 AliTPCFastVector & xxx = *xxxx;
2878 Int_t maxhits = kcmaxhits;
2881 AliTPChit * tpcHit=0;
2882 tpcHit = (AliTPChit*)FirstHit2(-1);
2883 Int_t currentIndex=0;
2884 Int_t lastrow=-1; //last writen row
2887 if (tpcHit==0) continue;
2888 sector=tpcHit->fSector; // sector number
2889 ipart=tpcHit->Track();
2893 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2894 Int_t index[3]={1,sector,0};
2895 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2896 if (currentrow<0) continue;
2897 if (lastrow<0) lastrow=currentrow;
2898 if (currentrow==lastrow){
2899 if ( currentIndex>=maxhits){
2901 xxx.ResizeTo(4*maxhits);
2903 xxx(currentIndex*4)=x[0];
2904 xxx(currentIndex*4+1)=x[1];
2905 xxx(currentIndex*4+2)=x[2];
2906 xxx(currentIndex*4+3)=tpcHit->fQ;
2910 if (currentIndex>2){
2922 for (Int_t index=0;index<currentIndex;index++){
2924 x=x2=x3=x4=xxx(index*4);
2932 sumy+=xxx(index*4+1);
2933 sumxy+=xxx(index*4+1)*x;
2934 sumx2y+=xxx(index*4+1)*x2;
2935 sumz+=xxx(index*4+2);
2936 sumxz+=xxx(index*4+2)*x;
2937 sumx2z+=xxx(index*4+2)*x2;
2938 sumq+=xxx(index*4+3);
2940 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
2941 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
2942 sumx2*(sumx*sumx3-sumx2*sumx2);
2944 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
2945 sumx2*(sumxy*sumx3-sumx2y*sumx2);
2946 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
2947 sumx2*(sumxz*sumx3-sumx2z*sumx2);
2949 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
2950 sumx2*(sumx*sumx2y-sumx2*sumxy);
2951 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
2952 sumx2*(sumx*sumx2z-sumx2*sumxz);
2954 Float_t y=detay/det+centralPad;
2955 Float_t z=detaz/det;
2956 Float_t by=detby/det; //y angle
2957 Float_t bz=detbz/det; //z angle
2958 sumy/=Float_t(currentIndex);
2959 sumz/=Float_t(currentIndex);
2961 AliComplexCluster cl;
2967 cl.fTracks[0]=ipart;
2969 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
2970 if (row!=0) row->InsertCluster(&cl);
2973 } //end of calculating cluster for given row
2975 } // end of loop over hits
2979 //_______________________________________________________________________________
2980 void AliTPC::Digits2Reco(Int_t firstevent,Int_t lastevent)
2982 // produces rec points from digits and writes them on the root file
2983 // AliTPCclusters.root
2985 TDirectory *cwd = gDirectory;
2988 AliTPCParamSR *dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60");
2990 printf("You are running 2 pad-length geom hits with 3 pad-length geom digits\n");
2992 dig = new AliTPCParamSR();
2996 dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");
2999 printf("No TPC parameters found\n");
3004 cout<<"AliTPC::Digits2Reco: TPC parameteres have been set"<<endl;
3006 if(!gSystem->Getenv("CONFIG_FILE")){
3007 out=TFile::Open("AliTPCclusters.root","recreate");
3013 tmp1=gSystem->Getenv("CONFIG_FILE_PREFIX");
3014 tmp2=gSystem->Getenv("CONFIG_OUTDIR");
3015 sprintf(tmp3,"%s%s/AliTPCclusters.root",tmp1,tmp2);
3016 out=TFile::Open(tmp3,"recreate");
3020 cout<<"AliTPC::Digits2Reco - determination of rec points begins"<<endl;
3023 for(Int_t iev=firstevent;iev<lastevent+1;iev++){
3025 printf("Processing event %d\n",iev);
3026 Digits2Clusters(iev);
3028 cout<<"AliTPC::Digits2Reco - determination of rec points ended"<<endl;
3037 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
3040 fLoader = new AliTPCLoader(GetName(),topfoldername);
3044 ////////////////////////////////////////////////////////////////////////
3045 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
3047 // load TPC paarmeters from a given file or create new if the object
3048 // is not found there
3049 // 12/05/2003 This method should be moved to the AliTPCLoader
3050 // and one has to decide where to store the TPC parameters
3053 sprintf(paramName,"75x40_100x60_150x60");
3054 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
3056 cout<<"TPC parameters "<<paramName<<" found."<<endl;
3058 cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3060 paramTPC = new AliTPCParamSR;
3064 // the older version of parameters can be accessed with this code.
3065 // In some cases, we have old parameters saved in the file but
3066 // digits were created with new parameters, it can be distinguish
3067 // by the name of TPC TreeD. The code here is just for the case
3068 // we would need to compare with old data, uncomment it if needed.
3070 // char paramName[50];
3071 // sprintf(paramName,"75x40_100x60");
3072 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
3074 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
3076 // sprintf(paramName,"75x40_100x60_150x60");
3077 // paramTPC=(AliTPCParam*)in->Get(paramName);
3079 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
3081 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3083 // paramTPC = new AliTPCParamSR;
3090 //_____________________________________________________________________________
3091 Double_t SigmaY2(Double_t r, Double_t tgl, Double_t pt)
3094 // Parametrised error of the cluster reconstruction (pad direction)
3097 const Float_t kArphi=0.41818e-2;
3098 const Float_t kBrphi=0.17460e-4;
3099 const Float_t kCrphi=0.30993e-2;
3100 const Float_t kDrphi=0.41061e-3;
3102 pt=TMath::Abs(pt)*1000.;
3104 tgl=TMath::Abs(tgl);
3105 Double_t s=kArphi - kBrphi*r*tgl + kCrphi*x*x + kDrphi*x;
3106 if (s<0.4e-3) s=0.4e-3;
3107 s*=1.3; //Iouri Belikov
3112 //_____________________________________________________________________________
3113 Double_t SigmaZ2(Double_t r, Double_t tgl)
3116 // Parametrised error of the cluster reconstruction (drift direction)
3119 const Float_t kAz=0.39614e-2;
3120 const Float_t kBz=0.22443e-4;
3121 const Float_t kCz=0.51504e-1;
3124 tgl=TMath::Abs(tgl);
3125 Double_t s=kAz - kBz*r*tgl + kCz*tgl*tgl;
3126 if (s<0.4e-3) s=0.4e-3;
3127 s*=1.3; //Iouri Belikov