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"
81 #include "AliTPCclusterer.h"
82 #include "AliTPCtracker.h"
83 #include "AliTPCpidESD.h"
88 //_____________________________________________________________________________
89 // helper class for fast matrix and vector manipulation - no range checking
90 // origin - Marian Ivanov
92 class AliTPCFastMatrix : public TMatrix {
94 AliTPCFastMatrix(Int_t rowlwb, Int_t rowupb, Int_t collwb, Int_t colupb);
95 Float_t & UncheckedAt(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
96 Float_t UncheckedAtFast(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
99 AliTPCFastMatrix::AliTPCFastMatrix(Int_t rowlwb, Int_t rowupb, Int_t collwb, Int_t colupb):
100 TMatrix(rowlwb, rowupb,collwb,colupb)
103 //_____________________________________________________________________________
104 class AliTPCFastVector : public TVector {
106 AliTPCFastVector(Int_t size);
107 virtual ~AliTPCFastVector(){;}
108 Float_t & UncheckedAt(Int_t index) const {return fElements[index];} //fast acces
112 AliTPCFastVector::AliTPCFastVector(Int_t size):TVector(size){
115 //_____________________________________________________________________________
119 // Default constructor
130 fHitType = 2; //default CONTAINERS - based on ROOT structure
137 //_____________________________________________________________________________
138 AliTPC::AliTPC(const char *name, const char *title)
139 : AliDetector(name,title)
142 // Standard constructor
146 // Initialise arrays of hits and digits
147 fHits = new TClonesArray("AliTPChit", 176);
148 gAlice->GetMCApp()->AddHitList(fHits);
153 fTrackHits = new AliTPCTrackHitsV2;
154 fTrackHits->SetHitPrecision(0.002);
155 fTrackHits->SetStepPrecision(0.003);
156 fTrackHits->SetMaxDistance(100);
158 fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
159 fTrackHitsOld->SetHitPrecision(0.002);
160 fTrackHitsOld->SetStepPrecision(0.003);
161 fTrackHitsOld->SetMaxDistance(100);
168 // Initialise counters
174 // Initialise color attributes
175 SetMarkerColor(kYellow);
178 // Set TPC parameters
182 if (!strcmp(title,"Default")) {
183 fTPCParam = new AliTPCParamSR;
185 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
191 //_____________________________________________________________________________
192 AliTPC::AliTPC(const AliTPC& t):AliDetector(t){
194 // dummy copy constructor
207 delete fTrackHits; //MI 15.09.2000
208 delete fTrackHitsOld; //MI 10.12.2001
209 if (fNoiseTable) delete [] fNoiseTable;
213 //_____________________________________________________________________________
214 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
217 // Add a hit to the list
219 // TClonesArray &lhits = *fHits;
220 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
222 TClonesArray &lhits = *fHits;
223 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
226 AddHit2(track,vol,hits);
229 //_____________________________________________________________________________
230 void AliTPC::BuildGeometry()
234 // Build TPC ROOT TNode geometry for the event display
239 const int kColorTPC=19;
240 char name[5], title[25];
241 const Double_t kDegrad=TMath::Pi()/180;
242 const Double_t kRaddeg=180./TMath::Pi();
245 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
246 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
248 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
249 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
251 Int_t nLo = fTPCParam->GetNInnerSector()/2;
252 Int_t nHi = fTPCParam->GetNOuterSector()/2;
254 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
255 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
256 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
257 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
260 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
261 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
267 // Get ALICE top node
270 nTop=gAlice->GetGeometry()->GetNode("alice");
274 rl = fTPCParam->GetInnerRadiusLow();
275 ru = fTPCParam->GetInnerRadiusUp();
279 sprintf(name,"LS%2.2d",i);
281 sprintf(title,"TPC low sector %3d",i);
284 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
285 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
286 tubs->SetNumberOfDivisions(1);
288 nNode = new TNode(name,title,name,0,0,0,"");
289 nNode->SetLineColor(kColorTPC);
295 rl = fTPCParam->GetOuterRadiusLow();
296 ru = fTPCParam->GetOuterRadiusUp();
299 sprintf(name,"US%2.2d",i);
301 sprintf(title,"TPC upper sector %d",i);
303 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
304 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
305 tubs->SetNumberOfDivisions(1);
307 nNode = new TNode(name,title,name,0,0,0,"");
308 nNode->SetLineColor(kColorTPC);
314 //_____________________________________________________________________________
315 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t ) const
318 // Calculate distance from TPC to mouse on the display
324 void AliTPC::Clusters2Tracks() const
326 //-----------------------------------------------------------------
327 // This is a track finder.
328 //-----------------------------------------------------------------
329 Error("Clusters2Tracks",
330 "Dummy function ! Call AliTPCtracker::Clusters2Tracks(...) instead !");
334 //_____________________________________________________________________________
335 void AliTPC::Reconstruct() const
337 // reconstruct clusters
339 AliLoader* loader = GetLoader();
340 loader->LoadRecPoints("recreate");
341 loader->LoadDigits("read");
343 AliTPCclusterer clusterer(fTPCParam);
344 AliRunLoader* runLoader = loader->GetRunLoader();
345 Int_t nEvents = runLoader->GetNumberOfEvents();
347 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
348 runLoader->GetEvent(iEvent);
350 TTree* treeClusters = loader->TreeR();
352 loader->MakeTree("R");
353 treeClusters = loader->TreeR();
355 TTree* treeDigits = loader->TreeD();
357 Error("Reconstruct", "Can't get digits tree !");
361 clusterer.Digits2Clusters(treeDigits, treeClusters);
363 loader->WriteRecPoints("OVERWRITE");
366 loader->UnloadRecPoints();
367 loader->UnloadDigits();
370 //_____________________________________________________________________________
371 AliTracker* AliTPC::CreateTracker() const
373 // create a TPC tracker
375 return new AliTPCtracker(fTPCParam);
378 //_____________________________________________________________________________
379 void AliTPC::FillESD(AliESD* esd) const
383 Double_t parTPC[] = {47., 0.10, 10.};
384 AliTPCpidESD tpcPID(parTPC);
389 //_____________________________________________________________________________
390 void AliTPC::CreateMaterials()
392 //-----------------------------------------------
393 // Create Materials for for TPC simulations
394 //-----------------------------------------------
396 //-----------------------------------------------------------------
397 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
398 //-----------------------------------------------------------------
400 Int_t iSXFLD=gAlice->Field()->Integ();
401 Float_t sXMGMX=gAlice->Field()->Max();
403 Float_t amat[5]; // atomic numbers
404 Float_t zmat[5]; // z
405 Float_t wmat[5]; // proportions
411 //***************** Gases *************************
413 //-------------------------------------------------
415 //-------------------------------------------------
426 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
436 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
439 //--------------------------------------------------------------
441 //--------------------------------------------------------------
458 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
460 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
475 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
477 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
493 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
495 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
497 //----------------------------------------------------------------
498 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
499 //----------------------------------------------------------------
505 Float_t rho,absl,x0,buf[1];
509 for(nc = 0;nc<fNoComp;nc++)
512 // retrive material constants
514 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,x0,absl,buf,nbuf);
519 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
521 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
522 density += fMixtProp[nc]*rho; // density of the mixture
526 // mixture proportions by weight!
528 for(nc = 0;nc<fNoComp;nc++)
531 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
533 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
534 apure[nnc] : amol[nnc])/am;
538 // Drift gases 1 - nonsensitive, 2 - sensitive
540 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
541 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
550 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
553 //----------------------------------------------------------------------
555 //----------------------------------------------------------------------
577 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
599 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
617 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
635 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
637 // SiO2 - used later for the glass fiber
649 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
658 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
667 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
676 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
694 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
713 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
715 // Epoxy - C14 H20 O3
732 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
740 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
744 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,x0,absl,buf,nbuf);
748 wmat[0]=0.644; // by weight!
751 density=0.5*(1.25+2.265);
753 AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
757 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,x0,absl,buf,nbuf);
759 wmat[0]=0.725; // by weight!
764 AliMixture(39,"G10",amat,zmat,density,2,wmat);
769 //----------------------------------------------------------
770 // tracking media for gases
771 //----------------------------------------------------------
773 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
774 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
775 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
776 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
778 //-----------------------------------------------------------
779 // tracking media for solids
780 //-----------------------------------------------------------
782 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
783 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
784 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
785 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
786 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
787 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
788 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
789 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
790 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
791 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
792 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
793 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
797 void AliTPC::GenerNoise(Int_t tablesize)
800 //Generate table with noise
807 if (fNoiseTable) delete[] fNoiseTable;
808 fNoiseTable = new Float_t[tablesize];
809 fNoiseDepth = tablesize;
810 fCurrentNoise =0; //!index of the noise in the noise table
812 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
813 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
816 Float_t AliTPC::GetNoise()
818 // get noise from table
819 // if ((fCurrentNoise%10)==0)
820 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
821 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
822 return fNoiseTable[fCurrentNoise++];
823 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
827 Bool_t AliTPC::IsSectorActive(Int_t sec) const
830 // check if the sector is active
831 if (!fActiveSectors) return kTRUE;
832 else return fActiveSectors[sec];
835 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
837 // activate interesting sectors
838 SetTreeAddress();//just for security
839 if (fActiveSectors) delete [] fActiveSectors;
840 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
841 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
842 for (Int_t i=0;i<n;i++)
843 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
847 void AliTPC::SetActiveSectors(Int_t flag)
850 // activate sectors which were hitted by tracks
852 SetTreeAddress();//just for security
853 if (fHitType==0) return; // if Clones hit - not short volume ID information
854 if (fActiveSectors) delete [] fActiveSectors;
855 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
857 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kTRUE;
860 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
864 Fatal("SetActiveSectors","Can not find TreeH in folder");
867 if (fHitType>1) branch = TreeH()->GetBranch("TPC2");
868 else branch = TreeH()->GetBranch("TPC");
869 Stat_t ntracks = TreeH()->GetEntries();
870 // loop over all hits
871 cout<<"\nAliTPC::SetActiveSectors(): Got "<<ntracks<<" tracks\n";
873 for(Int_t track=0;track<ntracks;track++)
877 if (fTrackHits && fHitType&4) {
878 TBranch * br1 = TreeH()->GetBranch("fVolumes");
879 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
880 br1->GetEvent(track);
881 br2->GetEvent(track);
882 Int_t *volumes = fTrackHits->GetVolumes();
883 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
884 fActiveSectors[volumes[j]]=kTRUE;
888 if (fTrackHitsOld && fHitType&2) {
889 TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
891 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
892 for (UInt_t j=0;j<ar->GetSize();j++){
893 fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
902 void AliTPC::Digits2Clusters(Int_t /*eventnumber*/) const
904 //-----------------------------------------------------------------
905 // This is a simple cluster finder.
906 //-----------------------------------------------------------------
907 Error("Digits2Clusters",
908 "Dummy function ! Call AliTPCclusterer::Digits2Clusters(...) instead !");
911 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
912 extern Double_t SigmaZ2(Double_t, Double_t);
913 //_____________________________________________________________________________
914 void AliTPC::Hits2Clusters(Int_t /*eventn*/)
916 //--------------------------------------------------------
917 // TPC simple cluster generator from hits
918 // obtained from the TPC Fast Simulator
919 // The point errors are taken from the parametrization
920 //--------------------------------------------------------
922 //-----------------------------------------------------------------
923 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
924 //-----------------------------------------------------------------
925 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
926 // Jouri.Belikov@cern.ch
927 //----------------------------------------------------------------
929 /////////////////////////////////////////////////////////////////////////////
931 //---------------------------------------------------------------------
932 // ALICE TPC Cluster Parameters
933 //--------------------------------------------------------------------
937 // Cluster width in rphi
938 const Float_t kACrphi=0.18322;
939 const Float_t kBCrphi=0.59551e-3;
940 const Float_t kCCrphi=0.60952e-1;
941 // Cluster width in z
942 const Float_t kACz=0.19081;
943 const Float_t kBCz=0.55938e-3;
944 const Float_t kCCz=0.30428;
948 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
952 //if(fDefaults == 0) SetDefaults();
954 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
956 TParticle *particle; // pointer to a given particle
957 AliTPChit *tpcHit; // pointer to a sigle TPC hit
961 Float_t pl,pt,tanth,rpad,ratio;
964 //---------------------------------------------------------------
965 // Get the access to the tracks
966 //---------------------------------------------------------------
971 Fatal("Hits2Clusters","Can not find TreeH in folder");
976 Stat_t ntracks = tH->GetEntries();
978 //Switch to the output file
980 if (fLoader->TreeR() == 0x0) fLoader->MakeTree("R");
982 cout<<"fTPCParam->GetTitle() = "<<fTPCParam->GetTitle()<<endl;
984 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
986 //fTPCParam->Write(fTPCParam->GetTitle());
988 AliTPCClustersArray carray;
989 carray.Setup(fTPCParam);
990 carray.SetClusterType("AliTPCcluster");
991 carray.MakeTree(fLoader->TreeR());
993 Int_t nclusters=0; //cluster counter
995 //------------------------------------------------------------
996 // Loop over all sectors (72 sectors for 20 deg
997 // segmentation for both lower and upper sectors)
998 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
999 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1001 // First cluster for sector 0 starts at "0"
1002 //------------------------------------------------------------
1004 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
1006 fTPCParam->AdjustCosSin(isec,cph,sph);
1008 //------------------------------------------------------------
1010 //------------------------------------------------------------
1012 for(Int_t track=0;track<ntracks;track++){
1014 tH->GetEvent(track);
1016 // Get number of the TPC hits
1018 tpcHit = (AliTPChit*)FirstHit(-1);
1024 if (tpcHit->fQ == 0.) {
1025 tpcHit = (AliTPChit*) NextHit();
1026 continue; //information about track (I.Belikov)
1028 sector=tpcHit->fSector; // sector number
1031 tpcHit = (AliTPChit*) NextHit();
1034 ipart=tpcHit->Track();
1035 particle=gAlice->GetMCApp()->Particle(ipart);
1038 if(pt < 1.e-9) pt=1.e-9;
1040 tanth = TMath::Abs(tanth);
1041 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
1042 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
1044 // space-point resolutions
1046 sigmaRphi=SigmaY2(rpad,tanth,pt);
1047 sigmaZ =SigmaZ2(rpad,tanth );
1051 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
1052 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
1054 // temporary protection
1056 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
1057 if(sigmaZ < 0.) sigmaZ=0.4e-3;
1058 if(clRphi < 0.) clRphi=2.5e-3;
1059 if(clZ < 0.) clZ=2.5e-5;
1064 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
1065 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1067 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
1068 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
1069 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
1070 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
1071 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1072 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
1073 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
1074 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
1075 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
1076 xyz[2]=sigmaRphi; // fSigmaY2
1077 xyz[3]=sigmaZ; // fSigmaZ2
1078 xyz[4]=tpcHit->fQ; // q
1080 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
1081 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
1083 Int_t tracks[3]={tpcHit->Track(), -1, -1};
1084 AliTPCcluster cluster(tracks,xyz);
1086 clrow->InsertCluster(&cluster); nclusters++;
1088 tpcHit = (AliTPChit*)NextHit();
1091 } // end of loop over hits
1093 } // end of loop over tracks
1095 Int_t nrows=fTPCParam->GetNRow(isec);
1096 for (Int_t irow=0; irow<nrows; irow++) {
1097 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
1098 if (!clrow) continue;
1099 carray.StoreRow(isec,irow);
1100 carray.ClearRow(isec,irow);
1103 } // end of loop over sectors
1105 cerr<<"Number of made clusters : "<<nclusters<<" \n";
1106 fLoader->WriteRecPoints("OVERWRITE");
1109 } // end of function
1111 //_________________________________________________________________
1112 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1114 //--------------------------------------------------------
1115 //calculate exact cross point of track and given pad row
1116 //resulting values are expressed in "digit" coordinata
1117 //--------------------------------------------------------
1119 //-----------------------------------------------------------------
1120 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1121 //-----------------------------------------------------------------
1123 if (fClustersArray==0){
1127 TParticle *particle; // pointer to a given particle
1128 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1129 // Int_t sector,nhits;
1131 const Int_t kcmaxhits=30000;
1132 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
1133 AliTPCFastVector & xxx = *xxxx;
1134 Int_t maxhits = kcmaxhits;
1135 //construct array for each padrow
1136 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
1137 fClustersArray->CreateRow(isec,i);
1139 //---------------------------------------------------------------
1140 // Get the access to the tracks
1141 //---------------------------------------------------------------
1143 TTree *tH = TreeH();
1146 Fatal("Hits2Clusters","Can not find TreeH in folder");
1151 Stat_t ntracks = tH->GetEntries();
1152 Int_t npart = gAlice->GetMCApp()->GetNtrack();
1155 if (fHitType>1) branch = tH->GetBranch("TPC2");
1156 else branch = tH->GetBranch("TPC");
1158 //------------------------------------------------------------
1160 //------------------------------------------------------------
1162 for(Int_t track=0;track<ntracks;track++){
1163 Bool_t isInSector=kTRUE;
1165 isInSector = TrackInVolume(isec,track);
1166 if (!isInSector) continue;
1168 branch->GetEntry(track); // get next track
1170 // Get number of the TPC hits and a pointer
1174 Int_t currentIndex=0;
1175 Int_t lastrow=-1; //last writen row
1179 tpcHit = (AliTPChit*)FirstHit(-1);
1182 Int_t sector=tpcHit->fSector; // sector number
1184 tpcHit = (AliTPChit*) NextHit();
1188 ipart=tpcHit->Track();
1189 if (ipart<npart) particle=gAlice->GetMCApp()->Particle(ipart);
1193 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1194 Int_t index[3]={1,isec,0};
1195 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1196 if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
1197 if (lastrow<0) lastrow=currentrow;
1198 if (currentrow==lastrow){
1199 if ( currentIndex>=maxhits){
1201 xxx.ResizeTo(4*maxhits);
1203 xxx(currentIndex*4)=x[0];
1204 xxx(currentIndex*4+1)=x[1];
1205 xxx(currentIndex*4+2)=x[2];
1206 xxx(currentIndex*4+3)=tpcHit->fQ;
1210 if (currentIndex>2){
1222 for (Int_t index=0;index<currentIndex;index++){
1224 x=x2=x3=x4=xxx(index*4);
1232 sumy+=xxx(index*4+1);
1233 sumxy+=xxx(index*4+1)*x;
1234 sumx2y+=xxx(index*4+1)*x2;
1235 sumz+=xxx(index*4+2);
1236 sumxz+=xxx(index*4+2)*x;
1237 sumx2z+=xxx(index*4+2)*x2;
1238 sumq+=xxx(index*4+3);
1240 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1241 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1242 sumx2*(sumx*sumx3-sumx2*sumx2);
1244 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1245 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1246 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1247 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1249 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1250 sumx2*(sumx*sumx2y-sumx2*sumxy);
1251 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1252 sumx2*(sumx*sumx2z-sumx2*sumxz);
1254 if (TMath::Abs(det)<0.00001){
1255 tpcHit = (AliTPChit*)NextHit();
1259 Float_t y=detay/det+centralPad;
1260 Float_t z=detaz/det;
1261 Float_t by=detby/det; //y angle
1262 Float_t bz=detbz/det; //z angle
1263 sumy/=Float_t(currentIndex);
1264 sumz/=Float_t(currentIndex);
1266 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1268 AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
1269 // AliComplexCluster cl;
1275 cl->fTracks[0]=ipart;
1279 } //end of calculating cluster for given row
1282 tpcHit = (AliTPChit*)NextHit();
1283 } // end of loop over hits
1284 } // end of loop over tracks
1285 //write padrows to tree
1286 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1287 fClustersArray->StoreRow(isec,ii);
1288 fClustersArray->ClearRow(isec,ii);
1296 //______________________________________________________________________
1297 AliDigitizer* AliTPC::CreateDigitizer(AliRunDigitizer* manager) const
1299 return new AliTPCDigitizer(manager);
1302 void AliTPC::SDigits2Digits2(Int_t /*eventnumber*/)
1304 //create digits from summable digits
1305 GenerNoise(500000); //create teble with noise
1307 //conect tree with sSDigits
1308 TTree *t = fLoader->TreeS();
1312 fLoader->LoadSDigits("READ");
1313 t = fLoader->TreeS();
1316 Error("SDigits2Digits2","Can not get input TreeS");
1321 if (fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
1323 AliSimDigits digarr, *dummy=&digarr;
1324 TBranch* sdb = t->GetBranch("Segment");
1327 Error("SDigits2Digits2","Can not find branch with segments in TreeS.");
1331 sdb->SetAddress(&dummy);
1333 Stat_t nentries = t->GetEntries();
1335 // set zero suppression
1337 fTPCParam->SetZeroSup(2);
1339 // get zero suppression
1341 Int_t zerosup = fTPCParam->GetZeroSup();
1343 //make tree with digits
1345 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1346 arr->SetClass("AliSimDigits");
1347 arr->Setup(fTPCParam);
1348 arr->MakeTree(fLoader->TreeD());
1350 AliTPCParam * par = fTPCParam;
1352 //Loop over segments of the TPC
1354 for (Int_t n=0; n<nentries; n++) {
1357 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1358 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
1361 if (!IsSectorActive(sec))
1363 // cout<<n<<" NOT Active \n";
1368 // cout<<n<<" Active \n";
1370 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1371 Int_t nrows = digrow->GetNRows();
1372 Int_t ncols = digrow->GetNCols();
1374 digrow->ExpandBuffer();
1375 digarr.ExpandBuffer();
1376 digrow->ExpandTrackBuffer();
1377 digarr.ExpandTrackBuffer();
1380 Short_t * pamp0 = digarr.GetDigits();
1381 Int_t * ptracks0 = digarr.GetTracks();
1382 Short_t * pamp1 = digrow->GetDigits();
1383 Int_t * ptracks1 = digrow->GetTracks();
1384 Int_t nelems =nrows*ncols;
1385 Int_t saturation = fTPCParam->GetADCSat();
1386 //use internal structure of the AliDigits - for speed reason
1387 //if you cahnge implementation
1388 //of the Alidigits - it must be rewriten -
1389 for (Int_t i= 0; i<nelems; i++){
1390 // Float_t q = *pamp0;
1391 //q/=16.; //conversion faktor
1392 //Float_t noise= GetNoise();
1394 //q= TMath::Nint(q);
1395 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1397 if (q>saturation) q=saturation;
1399 //if (ptracks0[0]==0)
1402 ptracks1[0]=ptracks0[0];
1403 ptracks1[nelems]=ptracks0[nelems];
1404 ptracks1[2*nelems]=ptracks0[2*nelems];
1412 arr->StoreRow(sec,row);
1413 arr->ClearRow(sec,row);
1414 // cerr<<sec<<"\t"<<row<<"\n";
1419 fLoader->WriteDigits("OVERWRITE");
1423 //__________________________________________________________________
1424 void AliTPC::SetDefaults(){
1426 // setting the defaults
1429 cerr<<"Setting default parameters...\n";
1431 // Set response functions
1434 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1436 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1438 printf("You are using 2 pad-length geom hits with 3 pad-lenght geom digits...\n");
1440 param = new AliTPCParamSR();
1443 param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
1446 printf("No TPC parameters found\n");
1451 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1452 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1453 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1454 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1455 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1456 rf->SetOffset(3*param->GetZSigma());
1459 TDirectory *savedir=gDirectory;
1460 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1462 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1467 prfinner->Read("prf_07504_Gati_056068_d02");
1468 //PH Set different names
1469 s=prfinner->GetGRF()->GetName();
1471 prfinner->GetGRF()->SetName(s.Data());
1473 prfouter1->Read("prf_10006_Gati_047051_d03");
1474 s=prfouter1->GetGRF()->GetName();
1476 prfouter1->GetGRF()->SetName(s.Data());
1478 prfouter2->Read("prf_15006_Gati_047051_d03");
1479 s=prfouter2->GetGRF()->GetName();
1481 prfouter2->GetGRF()->SetName(s.Data());
1486 param->SetInnerPRF(prfinner);
1487 param->SetOuter1PRF(prfouter1);
1488 param->SetOuter2PRF(prfouter2);
1489 param->SetTimeRF(rf);
1499 //__________________________________________________________________
1500 void AliTPC::Hits2Digits()
1503 // creates digits from hits
1506 fLoader->LoadHits("read");
1507 fLoader->LoadDigits("recreate");
1508 AliRunLoader* runLoader = fLoader->GetRunLoader();
1510 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1511 runLoader->GetEvent(iEvent);
1513 Hits2Digits(iEvent);
1516 fLoader->UnloadHits();
1517 fLoader->UnloadDigits();
1519 //__________________________________________________________________
1520 void AliTPC::Hits2Digits(Int_t eventnumber)
1522 //----------------------------------------------------
1523 // Loop over all sectors for a single event
1524 //----------------------------------------------------
1525 AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1526 rl->GetEvent(eventnumber);
1527 if (fLoader->TreeH() == 0x0)
1529 if(fLoader->LoadHits())
1531 Error("Hits2Digits","Can not load hits.");
1536 if (fLoader->TreeD() == 0x0 )
1538 fLoader->MakeTree("D");
1539 if (fLoader->TreeD() == 0x0 )
1541 Error("Hits2Digits","Can not get TreeD");
1546 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1547 GenerNoise(500000); //create teble with noise
1549 //setup TPCDigitsArray
1551 if(GetDigitsArray()) delete GetDigitsArray();
1553 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1554 arr->SetClass("AliSimDigits");
1555 arr->Setup(fTPCParam);
1557 arr->MakeTree(fLoader->TreeD());
1558 SetDigitsArray(arr);
1560 fDigitsSwitch=0; // standard digits
1562 cerr<<"Digitizing TPC -- normal digits...\n";
1564 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1565 if (IsSectorActive(isec))
1567 if (fDebug) Info("Hits2Digits","Sector %d is active.",isec);
1568 Hits2DigitsSector(isec);
1572 if (fDebug) Info("Hits2Digits","Sector %d is NOT active.",isec);
1575 fLoader->WriteDigits("OVERWRITE");
1577 //this line prevents the crash in the similar one
1578 //on the beginning of this method
1579 //destructor attempts to reset the tree, which is deleted by the loader
1580 //need to be redesign
1581 if(GetDigitsArray()) delete GetDigitsArray();
1582 SetDigitsArray(0x0);
1586 //__________________________________________________________________
1587 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1590 //-----------------------------------------------------------
1591 // summable digits - 16 bit "ADC", no noise, no saturation
1592 //-----------------------------------------------------------
1594 //----------------------------------------------------
1595 // Loop over all sectors for a single event
1596 //----------------------------------------------------
1597 // AliRunLoader* rl = (AliRunLoader*)fLoader->GetEventFolder()->FindObject(AliRunLoader::fgkRunLoaderName);
1599 AliRunLoader* rl = fLoader->GetRunLoader();
1601 rl->GetEvent(eventnumber);
1602 if (fLoader->TreeH() == 0x0)
1604 if(fLoader->LoadHits())
1606 Error("Hits2Digits","Can not load hits.");
1613 if (fLoader->TreeS() == 0x0 )
1615 fLoader->MakeTree("S");
1618 if(fDefaults == 0) SetDefaults();
1620 GenerNoise(500000); //create table with noise
1621 //setup TPCDigitsArray
1623 if(GetDigitsArray()) delete GetDigitsArray();
1626 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1627 arr->SetClass("AliSimDigits");
1628 arr->Setup(fTPCParam);
1629 arr->MakeTree(fLoader->TreeS());
1631 SetDigitsArray(arr);
1633 cerr<<"Digitizing TPC -- summable digits...\n";
1635 fDigitsSwitch=1; // summable digits
1637 // set zero suppression to "0"
1639 fTPCParam->SetZeroSup(0);
1641 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++)
1642 if (IsSectorActive(isec))
1644 // cout<<"Sector "<<isec<<" is active\n";
1645 Hits2DigitsSector(isec);
1648 fLoader->WriteSDigits("OVERWRITE");
1650 //this line prevents the crash in the similar one
1651 //on the beginning of this method
1652 //destructor attempts to reset the tree, which is deleted by the loader
1653 //need to be redesign
1654 if(GetDigitsArray()) delete GetDigitsArray();
1655 SetDigitsArray(0x0);
1657 //__________________________________________________________________
1659 void AliTPC::Hits2SDigits()
1662 //-----------------------------------------------------------
1663 // summable digits - 16 bit "ADC", no noise, no saturation
1664 //-----------------------------------------------------------
1666 fLoader->LoadHits("read");
1667 fLoader->LoadSDigits("recreate");
1668 AliRunLoader* runLoader = fLoader->GetRunLoader();
1670 for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
1671 runLoader->GetEvent(iEvent);
1674 Hits2SDigits2(iEvent);
1677 fLoader->UnloadHits();
1678 fLoader->UnloadSDigits();
1680 //_____________________________________________________________________________
1682 void AliTPC::Hits2DigitsSector(Int_t isec)
1684 //-------------------------------------------------------------------
1685 // TPC conversion from hits to digits.
1686 //-------------------------------------------------------------------
1688 //-----------------------------------------------------------------
1689 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1690 //-----------------------------------------------------------------
1692 //-------------------------------------------------------
1693 // Get the access to the track hits
1694 //-------------------------------------------------------
1696 // check if the parameters are set - important if one calls this method
1697 // directly, not from the Hits2Digits
1699 if(fDefaults == 0) SetDefaults();
1701 TTree *tH = TreeH(); // pointer to the hits tree
1704 Fatal("Hits2DigitsSector","Can not find TreeH in folder");
1708 Stat_t ntracks = tH->GetEntries();
1712 //-------------------------------------------
1713 // Only if there are any tracks...
1714 //-------------------------------------------
1718 //printf("*** Processing sector number %d ***\n",isec);
1720 Int_t nrows =fTPCParam->GetNRow(isec);
1722 row= new TObjArray* [nrows+2]; // 2 extra rows for cross talk
1724 MakeSector(isec,nrows,tH,ntracks,row);
1726 //--------------------------------------------------------
1727 // Digitize this sector, row by row
1728 // row[i] is the pointer to the TObjArray of AliTPCFastVectors,
1729 // each one containing electrons accepted on this
1730 // row, assigned into tracks
1731 //--------------------------------------------------------
1735 if (fDigitsArray->GetTree()==0)
1737 Fatal("Hits2DigitsSector","Tree not set in fDigitsArray");
1740 for (i=0;i<nrows;i++){
1742 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1744 DigitizeRow(i,isec,row);
1746 fDigitsArray->StoreRow(isec,i);
1748 Int_t ndig = dig->GetDigitSize();
1751 printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1753 fDigitsArray->ClearRow(isec,i);
1756 } // end of the sector digitization
1758 for(i=0;i<nrows+2;i++){
1763 delete [] row; // delete the array of pointers to TObjArray-s
1767 } // end of Hits2DigitsSector
1770 //_____________________________________________________________________________
1771 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1773 //-----------------------------------------------------------
1774 // Single row digitization, coupling from the neighbouring
1775 // rows taken into account
1776 //-----------------------------------------------------------
1778 //-----------------------------------------------------------------
1779 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1780 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1781 //-----------------------------------------------------------------
1784 Float_t zerosup = fTPCParam->GetZeroSup();
1785 // Int_t nrows =fTPCParam->GetNRow(isec);
1786 fCurrentIndex[1]= isec;
1789 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1790 Int_t nofTbins = fTPCParam->GetMaxTBin();
1791 Int_t indexRange[4];
1793 // Integrated signal for this row
1794 // and a single track signal
1797 AliTPCFastMatrix *m1 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // integrated
1798 AliTPCFastMatrix *m2 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // single
1800 AliTPCFastMatrix &total = *m1;
1802 // Array of pointers to the label-signal list
1804 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1805 Float_t **pList = new Float_t* [nofDigits];
1809 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1813 //Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1814 //Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1817 for (Int_t row= row1;row<=row2;row++){
1818 Int_t nTracks= rows[row]->GetEntries();
1819 for (i1=0;i1<nTracks;i1++){
1820 fCurrentIndex[2]= row;
1821 fCurrentIndex[3]=irow+1;
1823 m2->Zero(); // clear single track signal matrix
1824 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1825 GetList(trackLabel,nofPads,m2,indexRange,pList);
1827 else GetSignal(rows[row],i1,0,m1,indexRange);
1833 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1835 Float_t fzerosup = zerosup+0.5;
1836 for(Int_t it=0;it<nofTbins;it++){
1837 Float_t *pq = &(total.UncheckedAt(0,it));
1838 for(Int_t ip=0;ip<nofPads;ip++){
1842 if(fDigitsSwitch == 0){
1844 if(q <=fzerosup) continue; // do not fill zeros
1846 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1851 if(q <= 0.) continue; // do not fill zeros
1852 if(q>2000.) q=2000.;
1858 // "real" signal or electronic noise (list = -1)?
1861 for(Int_t j1=0;j1<3;j1++){
1862 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1867 <A NAME="AliDigits"></A>
1868 using of AliDigits object
1871 dig->SetDigitFast((Short_t)q,it,ip);
1872 if (fDigitsArray->IsSimulated())
1874 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1875 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1876 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1880 } // end of loop over time buckets
1881 } // end of lop over pads
1884 // This row has been digitized, delete nonused stuff
1887 for(lp=0;lp<nofDigits;lp++){
1888 if(pList[lp]) delete [] pList[lp];
1897 } // end of DigitizeRow
1899 //_____________________________________________________________________________
1901 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1902 AliTPCFastMatrix *m1, AliTPCFastMatrix *m2,Int_t *indexRange)
1905 //---------------------------------------------------------------
1906 // Calculates 2-D signal (pad,time) for a single track,
1907 // returns a pointer to the signal matrix and the track label
1908 // No digitization is performed at this level!!!
1909 //---------------------------------------------------------------
1911 //-----------------------------------------------------------------
1912 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1913 // Modified: Marian Ivanov
1914 //-----------------------------------------------------------------
1916 AliTPCFastVector *tv;
1918 tv = (AliTPCFastVector*)p1->At(ntr); // pointer to a track
1919 AliTPCFastVector &v = *tv;
1921 Float_t label = v(0);
1922 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3]-1)-1)/2;
1924 Int_t nElectrons = (tv->GetNrows()-1)/4;
1925 indexRange[0]=9999; // min pad
1926 indexRange[1]=-1; // max pad
1927 indexRange[2]=9999; //min time
1928 indexRange[3]=-1; // max time
1930 AliTPCFastMatrix &signal = *m1;
1931 AliTPCFastMatrix &total = *m2;
1933 // Loop over all electrons
1935 for(Int_t nel=0; nel<nElectrons; nel++){
1937 Float_t aval = v(idx+4);
1938 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1939 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1940 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
1942 Int_t *index = fTPCParam->GetResBin(0);
1943 Float_t *weight = & (fTPCParam->GetResWeight(0));
1945 if (n>0) for (Int_t i =0; i<n; i++){
1946 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
1949 Int_t time=index[2];
1950 Float_t qweight = *(weight)*eltoadcfac;
1952 if (m1!=0) signal.UncheckedAt(pad,time)+=qweight;
1953 total.UncheckedAt(pad,time)+=qweight;
1954 if (indexRange[0]>pad) indexRange[0]=pad;
1955 if (indexRange[1]<pad) indexRange[1]=pad;
1956 if (indexRange[2]>time) indexRange[2]=time;
1957 if (indexRange[3]<time) indexRange[3]=time;
1964 } // end of loop over electrons
1966 return label; // returns track label when finished
1969 //_____________________________________________________________________________
1970 void AliTPC::GetList(Float_t label,Int_t np,AliTPCFastMatrix *m,
1971 Int_t *indexRange, Float_t **pList)
1973 //----------------------------------------------------------------------
1974 // Updates the list of tracks contributing to digits for a given row
1975 //----------------------------------------------------------------------
1977 //-----------------------------------------------------------------
1978 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1979 //-----------------------------------------------------------------
1981 AliTPCFastMatrix &signal = *m;
1983 // lop over nonzero digits
1985 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
1986 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
1989 // accept only the contribution larger than 500 electrons (1/2 s_noise)
1991 if(signal(ip,it)<0.5) continue;
1994 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
1996 if(!pList[globalIndex]){
1999 // Create new list (6 elements - 3 signals and 3 labels),
2002 pList[globalIndex] = new Float_t [6];
2006 *pList[globalIndex] = -1.;
2007 *(pList[globalIndex]+1) = -1.;
2008 *(pList[globalIndex]+2) = -1.;
2009 *(pList[globalIndex]+3) = -1.;
2010 *(pList[globalIndex]+4) = -1.;
2011 *(pList[globalIndex]+5) = -1.;
2014 *pList[globalIndex] = label;
2015 *(pList[globalIndex]+3) = signal(ip,it);
2019 // check the signal magnitude
2021 Float_t highest = *(pList[globalIndex]+3);
2022 Float_t middle = *(pList[globalIndex]+4);
2023 Float_t lowest = *(pList[globalIndex]+5);
2026 // compare the new signal with already existing list
2029 if(signal(ip,it)<lowest) continue; // neglect this track
2033 if (signal(ip,it)>highest){
2034 *(pList[globalIndex]+5) = middle;
2035 *(pList[globalIndex]+4) = highest;
2036 *(pList[globalIndex]+3) = signal(ip,it);
2038 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2039 *(pList[globalIndex]+1) = *pList[globalIndex];
2040 *pList[globalIndex] = label;
2042 else if (signal(ip,it)>middle){
2043 *(pList[globalIndex]+5) = middle;
2044 *(pList[globalIndex]+4) = signal(ip,it);
2046 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2047 *(pList[globalIndex]+1) = label;
2050 *(pList[globalIndex]+5) = signal(ip,it);
2051 *(pList[globalIndex]+2) = label;
2055 } // end of loop over pads
2056 } // end of loop over time bins
2061 //___________________________________________________________________
2062 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
2063 Stat_t ntracks,TObjArray **row)
2066 //-----------------------------------------------------------------
2067 // Prepares the sector digitization, creates the vectors of
2068 // tracks for each row of this sector. The track vector
2069 // contains the track label and the position of electrons.
2070 //-----------------------------------------------------------------
2072 //-----------------------------------------------------------------
2073 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2074 //-----------------------------------------------------------------
2076 Float_t gasgain = fTPCParam->GetGasGain();
2080 AliTPChit *tpcHit; // pointer to a sigle TPC hit
2083 if (fHitType>1) branch = TH->GetBranch("TPC2");
2084 else branch = TH->GetBranch("TPC");
2087 //----------------------------------------------
2088 // Create TObjArray-s, one for each row,
2089 // each TObjArray will store the AliTPCFastVectors
2090 // of electrons, one AliTPCFastVectors per each track.
2091 //----------------------------------------------
2093 Int_t *nofElectrons = new Int_t [nrows+2]; // electron counter for each row
2094 AliTPCFastVector **tracks = new AliTPCFastVector* [nrows+2]; //pointers to the track vectors
2096 for(i=0; i<nrows+2; i++){
2097 row[i] = new TObjArray;
2104 //--------------------------------------------------------------------
2105 // Loop over tracks, the "track" contains the full history
2106 //--------------------------------------------------------------------
2108 Int_t previousTrack,currentTrack;
2109 previousTrack = -1; // nothing to store so far!
2111 for(Int_t track=0;track<ntracks;track++){
2112 Bool_t isInSector=kTRUE;
2114 isInSector = TrackInVolume(isec,track);
2115 if (!isInSector) continue;
2117 branch->GetEntry(track); // get next track
2121 tpcHit = (AliTPChit*)FirstHit(-1);
2123 //--------------------------------------------------------------
2125 //--------------------------------------------------------------
2130 Int_t sector=tpcHit->fSector; // sector number
2132 tpcHit = (AliTPChit*) NextHit();
2136 currentTrack = tpcHit->Track(); // track number
2139 if(currentTrack != previousTrack){
2141 // store already filled fTrack
2143 for(i=0;i<nrows+2;i++){
2144 if(previousTrack != -1){
2145 if(nofElectrons[i]>0){
2146 AliTPCFastVector &v = *tracks[i];
2147 v(0) = previousTrack;
2148 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2149 row[i]->Add(tracks[i]);
2152 delete tracks[i]; // delete empty AliTPCFastVector
2158 tracks[i] = new AliTPCFastVector(481); // AliTPCFastVectors for the next fTrack
2160 } // end of loop over rows
2162 previousTrack=currentTrack; // update track label
2165 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2167 //---------------------------------------------------
2168 // Calculate the electron attachment probability
2169 //---------------------------------------------------
2172 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
2173 /fTPCParam->GetDriftV();
2175 Float_t attProb = fTPCParam->GetAttCoef()*
2176 fTPCParam->GetOxyCont()*time; // fraction!
2178 //-----------------------------------------------
2179 // Loop over electrons
2180 //-----------------------------------------------
2183 for(Int_t nel=0;nel<qI;nel++){
2184 // skip if electron lost due to the attachment
2185 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2190 // protection for the nonphysical avalanche size (10**6 maximum)
2192 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2193 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2196 TransportElectron(xyz,index);
2198 fTPCParam->GetPadRow(xyz,index);
2199 // row 0 - cross talk from the innermost row
2200 // row fNRow+1 cross talk from the outermost row
2201 rowNumber = index[2]+1;
2202 //transform position to local digit coordinates
2203 //relative to nearest pad row
2204 if ((rowNumber<0)||rowNumber>fTPCParam->GetNRow(isec)+1) continue;
2206 if (isec <fTPCParam->GetNInnerSector()) {
2207 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2208 y1 = fTPCParam->GetYInner(rowNumber);
2211 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2212 y1 = fTPCParam->GetYOuter(rowNumber);
2214 // gain inefficiency at the wires edges - linear
2217 if(x1>y1) xyz[3]*=TMath::Max(1.e-6,(y1-x1+1.));
2219 nofElectrons[rowNumber]++;
2220 //----------------------------------
2221 // Expand vector if necessary
2222 //----------------------------------
2223 if(nofElectrons[rowNumber]>120){
2224 Int_t range = tracks[rowNumber]->GetNrows();
2225 if((nofElectrons[rowNumber])>(range-1)/4){
2227 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
2231 AliTPCFastVector &v = *tracks[rowNumber];
2232 Int_t idx = 4*nofElectrons[rowNumber]-3;
2233 Real_t * position = &(((AliTPCFastVector&)v).UncheckedAt(idx)); //make code faster
2234 memcpy(position,xyz,4*sizeof(Float_t));
2236 } // end of loop over electrons
2238 tpcHit = (AliTPChit*)NextHit();
2240 } // end of loop over hits
2241 } // end of loop over tracks
2244 // store remaining track (the last one) if not empty
2247 for(i=0;i<nrows+2;i++){
2248 if(nofElectrons[i]>0){
2249 AliTPCFastVector &v = *tracks[i];
2250 v(0) = previousTrack;
2251 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2252 row[i]->Add(tracks[i]);
2261 delete [] nofElectrons;
2264 } // end of MakeSector
2267 //_____________________________________________________________________________
2271 // Initialise TPC detector after definition of geometry
2276 printf("\n%s: ",ClassName());
2277 for(i=0;i<35;i++) printf("*");
2278 printf(" TPC_INIT ");
2279 for(i=0;i<35;i++) printf("*");
2280 printf("\n%s: ",ClassName());
2282 for(i=0;i<80;i++) printf("*");
2287 //_____________________________________________________________________________
2288 void AliTPC::MakeBranch(Option_t* option)
2291 // Create Tree branches for the TPC.
2293 if(GetDebug()) Info("MakeBranch","");
2294 Int_t buffersize = 4000;
2295 char branchname[10];
2296 sprintf(branchname,"%s",GetName());
2298 const char *h = strstr(option,"H");
2300 if ( h && (fHitType<=1) && (fHits == 0x0)) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2302 AliDetector::MakeBranch(option);
2304 const char *d = strstr(option,"D");
2306 if (fDigits && fLoader->TreeD() && d)
2308 MakeBranchInTree(gAlice->TreeD(), branchname, &fDigits, buffersize, 0);
2311 if (fHitType>1) MakeBranch2(option,0); // MI change 14.09.2000
2314 //_____________________________________________________________________________
2315 void AliTPC::ResetDigits()
2318 // Reset number of digits and the digits array for this detector
2321 if (fDigits) fDigits->Clear();
2324 //_____________________________________________________________________________
2325 void AliTPC::SetSecAL(Int_t sec)
2327 //---------------------------------------------------
2328 // Activate/deactivate selection for lower sectors
2329 //---------------------------------------------------
2331 //-----------------------------------------------------------------
2332 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2333 //-----------------------------------------------------------------
2337 //_____________________________________________________________________________
2338 void AliTPC::SetSecAU(Int_t sec)
2340 //----------------------------------------------------
2341 // Activate/deactivate selection for upper sectors
2342 //---------------------------------------------------
2344 //-----------------------------------------------------------------
2345 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2346 //-----------------------------------------------------------------
2350 //_____________________________________________________________________________
2351 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2353 //----------------------------------------
2354 // Select active lower sectors
2355 //----------------------------------------
2357 //-----------------------------------------------------------------
2358 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2359 //-----------------------------------------------------------------
2369 //_____________________________________________________________________________
2370 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2371 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2372 Int_t s11 , Int_t s12)
2374 //--------------------------------
2375 // Select active upper sectors
2376 //--------------------------------
2378 //-----------------------------------------------------------------
2379 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2380 //-----------------------------------------------------------------
2396 //_____________________________________________________________________________
2397 void AliTPC::SetSens(Int_t sens)
2400 //-------------------------------------------------------------
2401 // Activates/deactivates the sensitive strips at the center of
2402 // the pad row -- this is for the space-point resolution calculations
2403 //-------------------------------------------------------------
2405 //-----------------------------------------------------------------
2406 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2407 //-----------------------------------------------------------------
2413 void AliTPC::SetSide(Float_t side=0.)
2415 // choice of the TPC side
2420 //____________________________________________________________________________
2421 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2422 Float_t p2,Float_t p3)
2425 // gax mixture definition
2439 //_____________________________________________________________________________
2441 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2444 // electron transport taking into account:
2446 // 2.ExB at the wires
2447 // 3. nonisochronity
2449 // xyz and index must be already transformed to system 1
2452 fTPCParam->Transform1to2(xyz,index);
2455 Float_t driftl=xyz[2];
2456 if(driftl<0.01) driftl=0.01;
2457 driftl=TMath::Sqrt(driftl);
2458 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2459 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2460 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2461 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2462 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2466 if (fTPCParam->GetMWPCReadout()==kTRUE){
2467 Float_t dx = fTPCParam->Transform2to2NearestWire(xyz,index);
2468 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2470 //add nonisochronity (not implemented yet)
2473 ClassImp(AliTPCdigit)
2475 //_____________________________________________________________________________
2476 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2480 // Creates a TPC digit object
2482 fSector = digits[0];
2483 fPadRow = digits[1];
2486 fSignal = digits[4];
2492 //_____________________________________________________________________________
2493 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2497 // Creates a TPC hit object
2507 //________________________________________________________________________
2508 // Additional code because of the AliTPCTrackHitsV2
2510 void AliTPC::MakeBranch2(Option_t *option,const char */*file*/)
2513 // Create a new branch in the current Root Tree
2514 // The branch of fHits is automatically split
2515 // MI change 14.09.2000
2516 if(GetDebug()) Info("MakeBranch2","");
2517 if (fHitType<2) return;
2518 char branchname[10];
2519 sprintf(branchname,"%s2",GetName());
2521 // Get the pointer to the header
2522 const char *cH = strstr(option,"H");
2524 if (fTrackHits && TreeH() && cH && fHitType&4)
2526 if(GetDebug()) Info("MakeBranch2","Making branch for Type 4 Hits");
2527 TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,fBufferSize,99);
2530 if (fTrackHitsOld && TreeH() && cH && fHitType&2)
2532 if(GetDebug()) Info("MakeBranch2","Making branch for Type 2 Hits");
2533 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2534 TreeH(),fBufferSize,99);
2535 TreeH()->GetListOfBranches()->Add(branch);
2539 void AliTPC::SetTreeAddress()
2541 //Sets tree address for hits
2544 if (fHits == 0x0 ) fHits = new TClonesArray("AliTPChit", 176);//skowron 20.06.03
2545 AliDetector::SetTreeAddress();
2547 if (fHitType>1) SetTreeAddress2();
2550 void AliTPC::SetTreeAddress2()
2553 // Set branch address for the TrackHits Tree
2555 if(GetDebug()) Info("SetTreeAddress2","");
2558 char branchname[20];
2559 sprintf(branchname,"%s2",GetName());
2561 // Branch address for hit tree
2562 TTree *treeH = TreeH();
2563 if ((treeH)&&(fHitType&4)) {
2564 branch = treeH->GetBranch(branchname);
2567 branch->SetAddress(&fTrackHits);
2568 if (GetDebug()) Info("SetTreeAddress2","fHitType&4 Setting");
2571 if (GetDebug()) Info("SetTreeAddress2","fHitType&4 Failed (can not find branch)");
2574 if ((treeH)&&(fHitType&2)) {
2575 branch = treeH->GetBranch(branchname);
2578 branch->SetAddress(&fTrackHitsOld);
2579 if (GetDebug()) Info("SetTreeAddress2","fHitType&2 Setting");
2581 else if (GetDebug())
2582 Info("SetTreeAddress2","fHitType&2 Failed (can not find branch)");
2584 //set address to TREETR
2586 TTree *treeTR = TreeTR();
2587 if (treeTR && fTrackReferences) {
2588 branch = treeTR->GetBranch(GetName());
2589 if (branch) branch->SetAddress(&fTrackReferences);
2594 void AliTPC::FinishPrimary()
2596 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2597 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2601 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2604 // add hit to the list
2607 int primary = gAlice->GetMCApp()->GetPrimary(track);
2608 gAlice->GetMCApp()->Particle(primary)->SetBit(kKeepBit);
2612 gAlice->GetMCApp()->FlagTrack(track);
2614 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2615 //if (hit->fTrack!=rtrack)
2616 // cout<<"bad track number\n";
2617 if (fTrackHits && fHitType&4)
2618 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2619 hits[1],hits[2],(Int_t)hits[3]);
2620 if (fTrackHitsOld &&fHitType&2 )
2621 fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2622 hits[1],hits[2],(Int_t)hits[3]);
2626 void AliTPC::ResetHits()
2628 if (fHitType&1) AliDetector::ResetHits();
2629 if (fHitType>1) ResetHits2();
2632 void AliTPC::ResetHits2()
2636 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2637 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2641 AliHit* AliTPC::FirstHit(Int_t track)
2643 if (fHitType>1) return FirstHit2(track);
2644 return AliDetector::FirstHit(track);
2646 AliHit* AliTPC::NextHit()
2651 if (fHitType>1) return NextHit2();
2653 return AliDetector::NextHit();
2656 AliHit* AliTPC::FirstHit2(Int_t track)
2659 // Initialise the hit iterator
2660 // Return the address of the first hit for track
2661 // If track>=0 the track is read from disk
2662 // while if track<0 the first hit of the current
2663 // track is returned
2666 gAlice->ResetHits();
2667 TreeH()->GetEvent(track);
2670 if (fTrackHits && fHitType&4) {
2671 fTrackHits->First();
2672 return fTrackHits->GetHit();
2674 if (fTrackHitsOld && fHitType&2) {
2675 fTrackHitsOld->First();
2676 return fTrackHitsOld->GetHit();
2682 AliHit* AliTPC::NextHit2()
2685 //Return the next hit for the current track
2688 if (fTrackHitsOld && fHitType&2) {
2689 fTrackHitsOld->Next();
2690 return fTrackHitsOld->GetHit();
2694 return fTrackHits->GetHit();
2700 void AliTPC::LoadPoints(Int_t)
2704 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2707 if(fHitType==1) AliDetector::LoadPoints(a);
2708 else LoadPoints2(a);
2715 void AliTPC::RemapTrackHitIDs(Int_t *map)
2720 if (!fTrackHits) return;
2722 if (fTrackHitsOld && fHitType&2){
2723 AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2724 for (UInt_t i=0;i<arr->GetSize();i++){
2725 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2726 info->fTrackID = map[info->fTrackID];
2729 if (fTrackHitsOld && fHitType&4){
2730 TClonesArray * arr = fTrackHits->GetArray();;
2731 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2732 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2733 info->fTrackID = map[info->fTrackID];
2738 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2740 //return bool information - is track in given volume
2741 //load only part of the track information
2742 //return true if current track is in volume
2745 if (fTrackHitsOld && fHitType&2) {
2746 TBranch * br = TreeH()->GetBranch("fTrackHitsInfo");
2747 br->GetEvent(track);
2748 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2749 for (UInt_t j=0;j<ar->GetSize();j++){
2750 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2754 if (fTrackHits && fHitType&4) {
2755 TBranch * br1 = TreeH()->GetBranch("fVolumes");
2756 TBranch * br2 = TreeH()->GetBranch("fNVolumes");
2757 br2->GetEvent(track);
2758 br1->GetEvent(track);
2759 Int_t *volumes = fTrackHits->GetVolumes();
2760 Int_t nvolumes = fTrackHits->GetNVolumes();
2761 if (!volumes && nvolumes>0) {
2762 printf("Problematic track\t%d\t%d",track,nvolumes);
2765 for (Int_t j=0;j<nvolumes; j++)
2766 if (volumes[j]==id) return kTRUE;
2770 TBranch * br = TreeH()->GetBranch("fSector");
2771 br->GetEvent(track);
2772 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2773 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2780 //_____________________________________________________________________________
2781 void AliTPC::LoadPoints2(Int_t)
2784 // Store x, y, z of all hits in memory
2786 if (fTrackHits == 0 && fTrackHitsOld==0) return;
2789 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2790 if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2792 if (nhits == 0) return;
2793 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2794 if (fPoints == 0) fPoints = new TObjArray(tracks);
2797 Int_t *ntrk=new Int_t[tracks];
2798 Int_t *limi=new Int_t[tracks];
2799 Float_t **coor=new Float_t*[tracks];
2800 for(Int_t i=0;i<tracks;i++) {
2806 AliPoints *points = 0;
2809 Int_t chunk=nhits/4+1;
2811 // Loop over all the hits and store their position
2813 ahit = FirstHit2(-1);
2815 trk=ahit->GetTrack();
2816 if(ntrk[trk]==limi[trk]) {
2818 // Initialise a new track
2819 fp=new Float_t[3*(limi[trk]+chunk)];
2821 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2822 delete [] coor[trk];
2829 fp[3*ntrk[trk] ] = ahit->X();
2830 fp[3*ntrk[trk]+1] = ahit->Y();
2831 fp[3*ntrk[trk]+2] = ahit->Z();
2839 for(trk=0; trk<tracks; ++trk) {
2841 points = new AliPoints();
2842 points->SetMarkerColor(GetMarkerColor());
2843 points->SetMarkerSize(GetMarkerSize());
2844 points->SetDetector(this);
2845 points->SetParticle(trk);
2846 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2847 fPoints->AddAt(points,trk);
2848 delete [] coor[trk];
2858 //_____________________________________________________________________________
2859 void AliTPC::LoadPoints3(Int_t)
2862 // Store x, y, z of all hits in memory
2863 // - only intersection point with pad row
2864 if (fTrackHits == 0) return;
2866 Int_t nhits = fTrackHits->GetEntriesFast();
2867 if (nhits == 0) return;
2868 Int_t tracks = gAlice->GetMCApp()->GetNtrack();
2869 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2870 fPoints->Expand(2*tracks);
2873 Int_t *ntrk=new Int_t[tracks];
2874 Int_t *limi=new Int_t[tracks];
2875 Float_t **coor=new Float_t*[tracks];
2876 for(Int_t i=0;i<tracks;i++) {
2882 AliPoints *points = 0;
2885 Int_t chunk=nhits/4+1;
2887 // Loop over all the hits and store their position
2889 ahit = FirstHit2(-1);
2890 //for (Int_t hit=0;hit<nhits;hit++) {
2894 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2895 trk=ahit->GetTrack();
2896 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2897 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2898 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2899 if (currentrow!=lastrow){
2900 lastrow = currentrow;
2901 //later calculate intersection point
2902 if(ntrk[trk]==limi[trk]) {
2904 // Initialise a new track
2905 fp=new Float_t[3*(limi[trk]+chunk)];
2907 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2908 delete [] coor[trk];
2915 fp[3*ntrk[trk] ] = ahit->X();
2916 fp[3*ntrk[trk]+1] = ahit->Y();
2917 fp[3*ntrk[trk]+2] = ahit->Z();
2924 for(trk=0; trk<tracks; ++trk) {
2926 points = new AliPoints();
2927 points->SetMarkerColor(GetMarkerColor()+1);
2928 points->SetMarkerStyle(5);
2929 points->SetMarkerSize(0.2);
2930 points->SetDetector(this);
2931 points->SetParticle(trk);
2932 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2933 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2934 fPoints->AddAt(points,tracks+trk);
2935 delete [] coor[trk];
2946 void AliTPC::FindTrackHitsIntersection(TClonesArray * /*arr*/)
2950 //fill clones array with intersection of current point with the
2955 const Int_t kcmaxhits=30000;
2956 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
2957 AliTPCFastVector & xxx = *xxxx;
2958 Int_t maxhits = kcmaxhits;
2961 AliTPChit * tpcHit=0;
2962 tpcHit = (AliTPChit*)FirstHit2(-1);
2963 Int_t currentIndex=0;
2964 Int_t lastrow=-1; //last writen row
2967 if (tpcHit==0) continue;
2968 sector=tpcHit->fSector; // sector number
2969 ipart=tpcHit->Track();
2973 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
2974 Int_t index[3]={1,sector,0};
2975 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2976 if (currentrow<0) continue;
2977 if (lastrow<0) lastrow=currentrow;
2978 if (currentrow==lastrow){
2979 if ( currentIndex>=maxhits){
2981 xxx.ResizeTo(4*maxhits);
2983 xxx(currentIndex*4)=x[0];
2984 xxx(currentIndex*4+1)=x[1];
2985 xxx(currentIndex*4+2)=x[2];
2986 xxx(currentIndex*4+3)=tpcHit->fQ;
2990 if (currentIndex>2){
3002 for (Int_t index=0;index<currentIndex;index++){
3004 x=x2=x3=x4=xxx(index*4);
3012 sumy+=xxx(index*4+1);
3013 sumxy+=xxx(index*4+1)*x;
3014 sumx2y+=xxx(index*4+1)*x2;
3015 sumz+=xxx(index*4+2);
3016 sumxz+=xxx(index*4+2)*x;
3017 sumx2z+=xxx(index*4+2)*x2;
3018 sumq+=xxx(index*4+3);
3020 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
3021 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
3022 sumx2*(sumx*sumx3-sumx2*sumx2);
3024 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
3025 sumx2*(sumxy*sumx3-sumx2y*sumx2);
3026 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
3027 sumx2*(sumxz*sumx3-sumx2z*sumx2);
3029 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
3030 sumx2*(sumx*sumx2y-sumx2*sumxy);
3031 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
3032 sumx2*(sumx*sumx2z-sumx2*sumxz);
3034 Float_t y=detay/det+centralPad;
3035 Float_t z=detaz/det;
3036 Float_t by=detby/det; //y angle
3037 Float_t bz=detbz/det; //z angle
3038 sumy/=Float_t(currentIndex);
3039 sumz/=Float_t(currentIndex);
3041 AliComplexCluster cl;
3047 cl.fTracks[0]=ipart;
3049 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
3050 if (row!=0) row->InsertCluster(&cl);
3053 } //end of calculating cluster for given row
3055 } // end of loop over hits
3059 //_______________________________________________________________________________
3060 void AliTPC::Digits2Reco(Int_t firstevent,Int_t lastevent)
3062 // produces rec points from digits and writes them on the root file
3063 // AliTPCclusters.root
3065 TDirectory *cwd = gDirectory;
3068 AliTPCParamSR *dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60");
3070 printf("You are running 2 pad-length geom hits with 3 pad-length geom digits\n");
3072 dig = new AliTPCParamSR();
3076 dig=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");
3079 printf("No TPC parameters found\n");
3084 cout<<"AliTPC::Digits2Reco: TPC parameteres have been set"<<endl;
3086 if(!gSystem->Getenv("CONFIG_FILE")){
3087 out=TFile::Open("AliTPCclusters.root","recreate");
3093 tmp1=gSystem->Getenv("CONFIG_FILE_PREFIX");
3094 tmp2=gSystem->Getenv("CONFIG_OUTDIR");
3095 sprintf(tmp3,"%s%s/AliTPCclusters.root",tmp1,tmp2);
3096 out=TFile::Open(tmp3,"recreate");
3100 cout<<"AliTPC::Digits2Reco - determination of rec points begins"<<endl;
3103 for(Int_t iev=firstevent;iev<lastevent+1;iev++){
3105 printf("Processing event %d\n",iev);
3106 Digits2Clusters(iev);
3108 cout<<"AliTPC::Digits2Reco - determination of rec points ended"<<endl;
3117 AliLoader* AliTPC::MakeLoader(const char* topfoldername)
3120 fLoader = new AliTPCLoader(GetName(),topfoldername);
3124 ////////////////////////////////////////////////////////////////////////
3125 AliTPCParam* AliTPC::LoadTPCParam(TFile *file) {
3127 // load TPC paarmeters from a given file or create new if the object
3128 // is not found there
3129 // 12/05/2003 This method should be moved to the AliTPCLoader
3130 // and one has to decide where to store the TPC parameters
3133 sprintf(paramName,"75x40_100x60_150x60");
3134 AliTPCParam *paramTPC=(AliTPCParam*)file->Get(paramName);
3136 cout<<"TPC parameters "<<paramName<<" found."<<endl;
3138 cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3140 paramTPC = new AliTPCParamSR;
3144 // the older version of parameters can be accessed with this code.
3145 // In some cases, we have old parameters saved in the file but
3146 // digits were created with new parameters, it can be distinguish
3147 // by the name of TPC TreeD. The code here is just for the case
3148 // we would need to compare with old data, uncomment it if needed.
3150 // char paramName[50];
3151 // sprintf(paramName,"75x40_100x60");
3152 // AliTPCParam *paramTPC=(AliTPCParam*)in->Get(paramName);
3154 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
3156 // sprintf(paramName,"75x40_100x60_150x60");
3157 // paramTPC=(AliTPCParam*)in->Get(paramName);
3159 // cout<<"TPC parameters "<<paramName<<" found."<<endl;
3161 // cerr<<"TPC parameters not found. Create new (they may be incorrect)."
3163 // paramTPC = new AliTPCParamSR;
3171 //____________________________________________________________________________
3172 Double_t SigmaY2(Double_t r, Double_t tgl, Double_t pt)
3175 // Parametrised error of the cluster reconstruction (pad direction)
3178 const Float_t kArphi=0.41818e-2;
3179 const Float_t kBrphi=0.17460e-4;
3180 const Float_t kCrphi=0.30993e-2;
3181 const Float_t kDrphi=0.41061e-3;
3183 pt=TMath::Abs(pt)*1000.;
3185 tgl=TMath::Abs(tgl);
3186 Double_t s=kArphi - kBrphi*r*tgl + kCrphi*x*x + kDrphi*x;
3187 if (s<0.4e-3) s=0.4e-3;
3188 s*=1.3; //Iouri Belikov
3194 //____________________________________________________________________________
3195 Double_t SigmaZ2(Double_t r, Double_t tgl)
3198 // Parametrised error of the cluster reconstruction (drift direction)
3201 const Float_t kAz=0.39614e-2;
3202 const Float_t kBz=0.22443e-4;
3203 const Float_t kCz=0.51504e-1;
3206 tgl=TMath::Abs(tgl);
3207 Double_t s=kAz - kBz*r*tgl + kCz*tgl*tgl;
3208 if (s<0.4e-3) s=0.4e-3;
3209 s*=1.3; //Iouri Belikov