1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 Revision 1.53 2002/02/25 11:02:56 kowal2
19 Changes towards speeding up the code. Thanks to Marian Ivanov.
21 Revision 1.52 2002/02/18 09:26:09 kowal2
22 Removed compiler warning
24 Revision 1.51 2002/01/21 17:13:21 kowal2
25 New track hits using root containers. Setting active sectors added.
27 Revision 1.50 2001/12/06 14:16:19 kowal2
30 Revision 1.49 2001/11/30 11:55:37 hristov
31 Noise table created in Hits2SDigits (M.Ivanov)
33 Revision 1.48 2001/11/24 16:10:21 kowal2
36 Revision 1.47 2001/11/19 10:25:34 kowal2
37 Nearest integer instead of integer when converting to ADC counts
39 Revision 1.46 2001/11/07 06:47:12 kowal2
42 Revision 1.45 2001/11/03 13:33:48 kowal2
43 Updated algorithms in Hits2SDigits, SDigits2Digits,
47 Revision 1.44 2001/08/30 09:28:48 hristov
48 TTree names are explicitly set via SetName(name) and then Write() is called
50 Revision 1.43 2001/07/28 12:02:54 hristov
51 Branch split level set to 99
53 Revision 1.42 2001/07/28 11:38:52 hristov
54 Loop variable declared once
56 Revision 1.41 2001/07/28 10:53:50 hristov
57 Digitisation done according to the general scheme (M.Ivanov)
59 Revision 1.40 2001/07/27 13:03:14 hristov
60 Default Branch split level set to 99
62 Revision 1.39 2001/07/26 09:09:34 kowal2
63 Hits2Reco method added
65 Revision 1.38 2001/07/20 14:32:43 kowal2
66 Processing of many events possible now
68 Revision 1.37 2001/06/12 07:17:18 kowal2
69 Hits2SDigits method implemented (summable digits)
71 Revision 1.36 2001/05/16 14:57:25 alibrary
72 New files for folders and Stack
74 Revision 1.35 2001/05/08 16:02:22 kowal2
75 Updated material specifications
77 Revision 1.34 2001/05/08 15:00:15 hristov
78 Corrections for tracking in arbitrary magnenetic field. Changes towards a concept of global Alice track. Back propagation of reconstructed tracks (Yu.Belikov)
80 Revision 1.33 2001/04/03 12:40:43 kowal2
83 Revision 1.32 2001/03/12 17:47:36 hristov
84 Changes needed on Sun with CC 5.0
86 Revision 1.31 2001/03/12 08:21:50 kowal2
87 Corrected C++ bug in the material definitions
89 Revision 1.30 2001/03/01 17:34:47 kowal2
90 Correction due to the accuracy problem
92 Revision 1.29 2001/02/28 16:34:40 kowal2
93 Protection against nonphysical values of the avalanche size,
96 Revision 1.28 2001/01/26 19:57:19 hristov
97 Major upgrade of AliRoot code
99 Revision 1.27 2001/01/13 17:29:33 kowal2
100 Sun compiler correction
102 Revision 1.26 2001/01/10 07:59:43 kowal2
103 Corrections to load points from the noncompressed hits.
105 Revision 1.25 2000/11/02 07:25:31 kowal2
106 Changes due to the new hit structure.
109 Revision 1.24 2000/10/05 16:06:09 kowal2
110 Forward declarations. Changes due to a new class AliComplexCluster.
112 Revision 1.23 2000/10/02 21:28:18 fca
113 Removal of useless dependecies via forward declarations
115 Revision 1.22 2000/07/10 20:57:39 hristov
116 Update of TPC code and macros by M.Kowalski
118 Revision 1.19.2.4 2000/06/26 07:39:42 kowal2
119 Changes to obey the coding rules
121 Revision 1.19.2.3 2000/06/25 08:38:41 kowal2
122 Splitted from AliTPCtracking
124 Revision 1.19.2.2 2000/06/16 12:59:28 kowal2
125 Changed parameter settings
127 Revision 1.19.2.1 2000/06/09 07:15:07 kowal2
129 Defaults loaded automatically (hard-wired)
130 Optional parameters can be set via macro called in the constructor
132 Revision 1.19 2000/04/18 19:00:59 fca
133 Small bug fixes to TPC files
135 Revision 1.18 2000/04/17 09:37:33 kowal2
136 removed obsolete AliTPCDigitsDisplay.C
138 Revision 1.17.2.2 2000/04/10 08:15:12 kowal2
140 New, experimental data structure from M. Ivanov
141 New tracking algorithm
142 Different pad geometry for different sectors
143 Digitization rewritten
145 Revision 1.17.2.1 2000/04/10 07:56:53 kowal2
146 Not used anymore - removed
148 Revision 1.17 2000/01/19 17:17:30 fca
149 Introducing a list of lists of hits -- more hits allowed for detector now
151 Revision 1.16 1999/11/05 09:29:23 fca
152 Accept only signals > 0
154 Revision 1.15 1999/10/08 06:26:53 fca
155 Removed ClustersIndex - not used anymore
157 Revision 1.14 1999/09/29 09:24:33 fca
158 Introduction of the Copyright and cvs Log
162 ///////////////////////////////////////////////////////////////////////////////
164 // Time Projection Chamber //
165 // This class contains the basic functions for the Time Projection Chamber //
166 // detector. Functions specific to one particular geometry are //
167 // contained in the derived classes //
171 <img src="picts/AliTPCClass.gif">
176 ///////////////////////////////////////////////////////////////////////////////
184 #include <TGeometry.h>
187 #include <TObjectTable.h>
188 #include "TParticle.h"
194 #include <iostream.h>
201 #include "AliTPCParamSR.h"
202 #include "AliTPCPRF2D.h"
203 #include "AliTPCRF1D.h"
204 #include "AliDigits.h"
205 #include "AliSimDigits.h"
206 #include "AliTPCTrackHits.h"
207 #include "AliTPCTrackHitsV2.h"
208 #include "AliPoints.h"
209 #include "AliArrayBranch.h"
212 #include "AliTPCDigitsArray.h"
213 #include "AliComplexCluster.h"
214 #include "AliClusters.h"
215 #include "AliTPCClustersRow.h"
216 #include "AliTPCClustersArray.h"
218 #include "AliTPCcluster.h"
219 #include "AliTPCclusterer.h"
220 #include "AliTPCtracker.h"
222 #include <TInterpreter.h>
229 //_____________________________________________________________________________
230 // helper class for fast matrix and vector manipulation - no range checking
231 // origin - Marian Ivanov
233 class AliTPCFastMatrix : public TMatrix {
235 AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb);
236 inline Float_t & UncheckedAt(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
237 inline Float_t UncheckedAtFast(Int_t rown, Int_t coln) const {return (fIndex[coln])[rown];} //fast acces
240 AliTPCFastMatrix::AliTPCFastMatrix(Int_t row_lwb, Int_t row_upb, Int_t col_lwb, Int_t col_upb):
241 TMatrix(row_lwb, row_upb,col_lwb,col_upb)
244 //_____________________________________________________________________________
245 class AliTPCFastVector : public TVector {
247 AliTPCFastVector(Int_t size);
248 inline Float_t & UncheckedAt(Int_t index) const {return fElements[index];} //fast acces
251 AliTPCFastVector::AliTPCFastVector(Int_t size):TVector(size){
254 //_____________________________________________________________________________
258 // Default constructor
269 fHitType = 4; //default CONTAINERS - based on ROOT structure
276 //_____________________________________________________________________________
277 AliTPC::AliTPC(const char *name, const char *title)
278 : AliDetector(name,title)
281 // Standard constructor
285 // Initialise arrays of hits and digits
286 fHits = new TClonesArray("AliTPChit", 176);
287 gAlice->AddHitList(fHits);
292 fTrackHits = new AliTPCTrackHitsV2;
293 fTrackHits->SetHitPrecision(0.002);
294 fTrackHits->SetStepPrecision(0.003);
295 fTrackHits->SetMaxDistance(100);
297 fTrackHitsOld = new AliTPCTrackHits; //MI - 13.09.2000
298 fTrackHitsOld->SetHitPrecision(0.002);
299 fTrackHitsOld->SetStepPrecision(0.003);
300 fTrackHitsOld->SetMaxDistance(100);
307 // Initialise counters
313 // Initialise color attributes
314 SetMarkerColor(kYellow);
317 // Set TPC parameters
321 if (!strcmp(title,"Default")) {
322 fTPCParam = new AliTPCParamSR;
324 cerr<<"AliTPC warning: in Config.C you must set non-default parameters\n";
330 //_____________________________________________________________________________
341 delete fTrackHits; //MI 15.09.2000
342 delete fTrackHitsOld; //MI 10.12.2001
343 if (fNoiseTable) delete [] fNoiseTable;
347 //_____________________________________________________________________________
348 void AliTPC::AddHit(Int_t track, Int_t *vol, Float_t *hits)
351 // Add a hit to the list
353 // TClonesArray &lhits = *fHits;
354 // new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
356 TClonesArray &lhits = *fHits;
357 new(lhits[fNhits++]) AliTPChit(fIshunt,track,vol,hits);
360 AddHit2(track,vol,hits);
363 //_____________________________________________________________________________
364 void AliTPC::BuildGeometry()
368 // Build TPC ROOT TNode geometry for the event display
373 const int kColorTPC=19;
374 char name[5], title[25];
375 const Double_t kDegrad=TMath::Pi()/180;
376 const Double_t kRaddeg=180./TMath::Pi();
379 Float_t innerOpenAngle = fTPCParam->GetInnerAngle();
380 Float_t outerOpenAngle = fTPCParam->GetOuterAngle();
382 Float_t innerAngleShift = fTPCParam->GetInnerAngleShift();
383 Float_t outerAngleShift = fTPCParam->GetOuterAngleShift();
385 Int_t nLo = fTPCParam->GetNInnerSector()/2;
386 Int_t nHi = fTPCParam->GetNOuterSector()/2;
388 const Double_t kloAng = (Double_t)TMath::Nint(innerOpenAngle*kRaddeg);
389 const Double_t khiAng = (Double_t)TMath::Nint(outerOpenAngle*kRaddeg);
390 const Double_t kloAngSh = (Double_t)TMath::Nint(innerAngleShift*kRaddeg);
391 const Double_t khiAngSh = (Double_t)TMath::Nint(outerAngleShift*kRaddeg);
394 const Double_t kloCorr = 1/TMath::Cos(0.5*kloAng*kDegrad);
395 const Double_t khiCorr = 1/TMath::Cos(0.5*khiAng*kDegrad);
401 // Get ALICE top node
404 nTop=gAlice->GetGeometry()->GetNode("alice");
408 rl = fTPCParam->GetInnerRadiusLow();
409 ru = fTPCParam->GetInnerRadiusUp();
413 sprintf(name,"LS%2.2d",i);
415 sprintf(title,"TPC low sector %3d",i);
418 tubs = new TTUBS(name,title,"void",rl*kloCorr,ru*kloCorr,250.,
419 kloAng*(i-0.5)+kloAngSh,kloAng*(i+0.5)+kloAngSh);
420 tubs->SetNumberOfDivisions(1);
422 nNode = new TNode(name,title,name,0,0,0,"");
423 nNode->SetLineColor(kColorTPC);
429 rl = fTPCParam->GetOuterRadiusLow();
430 ru = fTPCParam->GetOuterRadiusUp();
433 sprintf(name,"US%2.2d",i);
435 sprintf(title,"TPC upper sector %d",i);
437 tubs = new TTUBS(name,title,"void",rl*khiCorr,ru*khiCorr,250,
438 khiAng*(i-0.5)+khiAngSh,khiAng*(i+0.5)+khiAngSh);
439 tubs->SetNumberOfDivisions(1);
441 nNode = new TNode(name,title,name,0,0,0,"");
442 nNode->SetLineColor(kColorTPC);
448 //_____________________________________________________________________________
449 Int_t AliTPC::DistancetoPrimitive(Int_t , Int_t )
452 // Calculate distance from TPC to mouse on the display
458 void AliTPC::Clusters2Tracks(TFile *of) {
459 //-----------------------------------------------------------------
460 // This is a track finder.
461 //-----------------------------------------------------------------
462 AliTPCtracker tracker(fTPCParam);
463 tracker.Clusters2Tracks(gFile,of);
466 //_____________________________________________________________________________
467 void AliTPC::CreateMaterials()
469 //-----------------------------------------------
470 // Create Materials for for TPC simulations
471 //-----------------------------------------------
473 //-----------------------------------------------------------------
474 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
475 //-----------------------------------------------------------------
477 Int_t iSXFLD=gAlice->Field()->Integ();
478 Float_t sXMGMX=gAlice->Field()->Max();
480 Float_t amat[5]; // atomic numbers
481 Float_t zmat[5]; // z
482 Float_t wmat[5]; // proportions
488 //***************** Gases *************************
490 //-------------------------------------------------
492 //-------------------------------------------------
503 AliMaterial(20,"Ne",amat[0],zmat[0],density,999.,999.);
513 AliMaterial(21,"Ar",amat[0],zmat[0],density,999.,999.);
516 //--------------------------------------------------------------
518 //--------------------------------------------------------------
535 amol[0] = amat[0]*wmat[0]+amat[1]*wmat[1];
537 AliMixture(10,"CO2",amat,zmat,density,-2,wmat);
552 amol[1] = amat[0]*wmat[0]+amat[1]*wmat[1];
554 AliMixture(11,"CF4",amat,zmat,density,-2,wmat);
570 amol[2] = amat[0]*wmat[0]+amat[1]*wmat[1];
572 AliMixture(12,"CH4",amat,zmat,density,-2,wmat);
574 //----------------------------------------------------------------
575 // gases - mixtures, ID >= 20 pure gases, <= 10 ID < 20 -compounds
576 //----------------------------------------------------------------
582 Float_t rho,absl,X0,buf[1];
586 for(nc = 0;nc<fNoComp;nc++)
589 // retrive material constants
591 gMC->Gfmate((*fIdmate)[fMixtComp[nc]],namate,a,z,rho,X0,absl,buf,nbuf);
596 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
598 am += fMixtProp[nc]*((fMixtComp[nc]>=20) ? apure[nnc] : amol[nnc]);
599 density += fMixtProp[nc]*rho; // density of the mixture
603 // mixture proportions by weight!
605 for(nc = 0;nc<fNoComp;nc++)
608 Int_t nnc = (fMixtComp[nc]>=20) ? fMixtComp[nc]%20 : fMixtComp[nc]%10;
610 wmat[nc] = fMixtProp[nc]*((fMixtComp[nc]>=20) ?
611 apure[nnc] : amol[nnc])/am;
615 // Drift gases 1 - nonsensitive, 2 - sensitive
617 AliMixture(31,"Drift gas 1",amat,zmat,density,fNoComp,wmat);
618 AliMixture(32,"Drift gas 2",amat,zmat,density,fNoComp,wmat);
627 AliMaterial(24,"Air",amat[0],zmat[0],density,999.,999.);
630 //----------------------------------------------------------------------
632 //----------------------------------------------------------------------
654 AliMixture(34,"Kevlar",amat,zmat,density,-4,wmat);
676 AliMixture(35,"NOMEX",amat,zmat,density,-4,wmat);
694 AliMixture(36,"Makrolon",amat,zmat,density,-3,wmat);
712 AliMixture(37, "Mylar",amat,zmat,density,-3,wmat);
714 // SiO2 - used later for the glass fiber
726 AliMixture(38,"SiO2",amat,zmat,2.2,-2,wmat); //SiO2 - quartz (rho=2.2)
735 AliMaterial(40,"Al",amat[0],zmat[0],density,999.,999.);
744 AliMaterial(41,"Si",amat[0],zmat[0],density,999.,999.);
753 AliMaterial(42,"Cu",amat[0],zmat[0],density,999.,999.);
771 AliMixture(43, "Tedlar",amat,zmat,density,-3,wmat);
790 AliMixture(44,"Plexiglas",amat,zmat,density,-3,wmat);
792 // Epoxy - C14 H20 O3
809 AliMixture(45,"Epoxy",amat,zmat,density,-3,wmat);
817 AliMaterial(46,"C",amat[0],zmat[0],density,999.,999.);
821 gMC->Gfmate((*fIdmate)[45],namate,amat[1],zmat[1],rho,X0,absl,buf,nbuf);
825 wmat[0]=0.644; // by weight!
828 density=0.5*(1.25+2.265);
830 AliMixture(47,"Cfiber",amat,zmat,density,2,wmat);
834 gMC->Gfmate((*fIdmate)[38],namate,amat[0],zmat[0],rho,X0,absl,buf,nbuf);
836 wmat[0]=0.725; // by weight!
841 AliMixture(39,"G10",amat,zmat,density,2,wmat);
846 //----------------------------------------------------------
847 // tracking media for gases
848 //----------------------------------------------------------
850 AliMedium(0, "Air", 24, 0, iSXFLD, sXMGMX, 10., 999., .1, .01, .1);
851 AliMedium(1, "Drift gas 1", 31, 0, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
852 AliMedium(2, "Drift gas 2", 32, 1, iSXFLD, sXMGMX, 10., 999.,.1,.001, .001);
853 AliMedium(3,"CO2",10,0, iSXFLD, sXMGMX, 10., 999.,.1, .001, .001);
855 //-----------------------------------------------------------
856 // tracking media for solids
857 //-----------------------------------------------------------
859 AliMedium(4,"Al",40,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
860 AliMedium(5,"Kevlar",34,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
861 AliMedium(6,"Nomex",35,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
862 AliMedium(7,"Makrolon",36,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
863 AliMedium(8,"Mylar",37,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
864 AliMedium(9,"Tedlar",43,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
865 AliMedium(10,"Cu",42,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
866 AliMedium(11,"Si",41,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
867 AliMedium(12,"G10",39,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
868 AliMedium(13,"Plexiglas",44,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
869 AliMedium(14,"Epoxy",45,0, iSXFLD, sXMGMX, 10., 999., .1, .0005, .001);
870 AliMedium(15,"Cfiber",47,0, iSXFLD, sXMGMX, 10., 999., .1, .001, .001);
874 void AliTPC::GenerNoise(Int_t tablesize)
877 //Generate table with noise
884 if (fNoiseTable) delete[] fNoiseTable;
885 fNoiseTable = new Float_t[tablesize];
886 fNoiseDepth = tablesize;
887 fCurrentNoise =0; //!index of the noise in the noise table
889 Float_t norm = fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac();
890 for (Int_t i=0;i<tablesize;i++) fNoiseTable[i]= gRandom->Gaus(0,norm);
893 Float_t AliTPC::GetNoise()
895 // get noise from table
896 // if ((fCurrentNoise%10)==0)
897 // fCurrentNoise= gRandom->Rndm()*fNoiseDepth;
898 if (fCurrentNoise>=fNoiseDepth) fCurrentNoise=0;
899 return fNoiseTable[fCurrentNoise++];
900 //gRandom->Gaus(0, fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac());
904 Bool_t AliTPC::IsSectorActive(Int_t sec)
907 // check if the sector is active
908 if (!fActiveSectors) return kTRUE;
909 else return fActiveSectors[sec];
912 void AliTPC::SetActiveSectors(Int_t * sectors, Int_t n)
914 // activate interesting sectors
915 if (fActiveSectors) delete [] fActiveSectors;
916 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
917 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
918 for (Int_t i=0;i<n;i++)
919 if ((sectors[i]>=0) && sectors[i]<fTPCParam->GetNSector()) fActiveSectors[sectors[i]]=kTRUE;
923 void AliTPC::SetActiveSectors()
926 // activate sectors which were hitted by tracks
928 if (fHitType==0) return; // if Clones hit - not short volume ID information
929 if (fActiveSectors) delete [] fActiveSectors;
930 fActiveSectors = new Bool_t[fTPCParam->GetNSector()];
931 for (Int_t i=0;i<fTPCParam->GetNSector();i++) fActiveSectors[i]=kFALSE;
933 if (fHitType>1) branch = gAlice->TreeH()->GetBranch("TPC2");
934 else branch = gAlice->TreeH()->GetBranch("TPC");
935 Stat_t ntracks = gAlice->TreeH()->GetEntries();
936 // loop over all hits
937 for(Int_t track=0;track<ntracks;track++){
940 if (fTrackHits && fHitType&4) {
941 TBranch * br1 = gAlice->TreeH()->GetBranch("fVolumes");
942 TBranch * br2 = gAlice->TreeH()->GetBranch("fNVolumes");
943 br1->GetEvent(track);
944 br2->GetEvent(track);
945 Int_t *volumes = fTrackHits->GetVolumes();
946 for (Int_t j=0;j<fTrackHits->GetNVolumes(); j++)
947 fActiveSectors[volumes[j]]=kTRUE;
951 if (fTrackHitsOld && fHitType&2) {
952 TBranch * br = gAlice->TreeH()->GetBranch("fTrackHitsInfo");
954 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
955 for (UInt_t j=0;j<ar->GetSize();j++){
956 fActiveSectors[((AliTrackHitsInfo*)ar->At(j))->fVolumeID] =kTRUE;
966 void AliTPC::Digits2Clusters(TFile *of, Int_t eventnumber)
968 //-----------------------------------------------------------------
969 // This is a simple cluster finder.
970 //-----------------------------------------------------------------
971 AliTPCclusterer::Digits2Clusters(fTPCParam,of,eventnumber);
974 extern Double_t SigmaY2(Double_t, Double_t, Double_t);
975 extern Double_t SigmaZ2(Double_t, Double_t);
976 //_____________________________________________________________________________
977 void AliTPC::Hits2Clusters(TFile *of, Int_t eventn)
979 //--------------------------------------------------------
980 // TPC simple cluster generator from hits
981 // obtained from the TPC Fast Simulator
982 // The point errors are taken from the parametrization
983 //--------------------------------------------------------
985 //-----------------------------------------------------------------
986 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
987 //-----------------------------------------------------------------
988 // Adopted to Marian's cluster data structure by I.Belikov, CERN,
989 // Jouri.Belikov@cern.ch
990 //----------------------------------------------------------------
992 /////////////////////////////////////////////////////////////////////////////
994 //---------------------------------------------------------------------
995 // ALICE TPC Cluster Parameters
996 //--------------------------------------------------------------------
1000 // Cluster width in rphi
1001 const Float_t kACrphi=0.18322;
1002 const Float_t kBCrphi=0.59551e-3;
1003 const Float_t kCCrphi=0.60952e-1;
1004 // Cluster width in z
1005 const Float_t kACz=0.19081;
1006 const Float_t kBCz=0.55938e-3;
1007 const Float_t kCCz=0.30428;
1009 TDirectory *savedir=gDirectory;
1011 if (!of->IsOpen()) {
1012 cerr<<"AliTPC::Hits2Clusters(): output file not open !\n";
1016 //if(fDefaults == 0) SetDefaults();
1018 Float_t sigmaRphi,sigmaZ,clRphi,clZ;
1020 TParticle *particle; // pointer to a given particle
1021 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1025 Float_t pl,pt,tanth,rpad,ratio;
1028 //---------------------------------------------------------------
1029 // Get the access to the tracks
1030 //---------------------------------------------------------------
1032 TTree *tH = gAlice->TreeH();
1033 Stat_t ntracks = tH->GetEntries();
1035 //Switch to the output file
1040 sprintf(cname,"TreeC_TPC_%d",eventn);
1042 fTPCParam->Write(fTPCParam->GetTitle());
1043 AliTPCClustersArray carray;
1044 carray.Setup(fTPCParam);
1045 carray.SetClusterType("AliTPCcluster");
1048 Int_t nclusters=0; //cluster counter
1050 //------------------------------------------------------------
1051 // Loop over all sectors (72 sectors for 20 deg
1052 // segmentation for both lower and upper sectors)
1053 // Sectors 0-35 are lower sectors, 0-17 z>0, 17-35 z<0
1054 // Sectors 36-71 are upper sectors, 36-53 z>0, 54-71 z<0
1056 // First cluster for sector 0 starts at "0"
1057 //------------------------------------------------------------
1059 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++){
1061 fTPCParam->AdjustCosSin(isec,cph,sph);
1063 //------------------------------------------------------------
1065 //------------------------------------------------------------
1067 for(Int_t track=0;track<ntracks;track++){
1069 tH->GetEvent(track);
1071 // Get number of the TPC hits
1073 tpcHit = (AliTPChit*)FirstHit(-1);
1079 if (tpcHit->fQ == 0.) {
1080 tpcHit = (AliTPChit*) NextHit();
1081 continue; //information about track (I.Belikov)
1083 sector=tpcHit->fSector; // sector number
1086 tpcHit = (AliTPChit*) NextHit();
1089 ipart=tpcHit->Track();
1090 particle=gAlice->Particle(ipart);
1093 if(pt < 1.e-9) pt=1.e-9;
1095 tanth = TMath::Abs(tanth);
1096 rpad=TMath::Sqrt(tpcHit->X()*tpcHit->X() + tpcHit->Y()*tpcHit->Y());
1097 ratio=0.001*rpad/pt; // pt must be in MeV/c - historical reason
1099 // space-point resolutions
1101 sigmaRphi=SigmaY2(rpad,tanth,pt);
1102 sigmaZ =SigmaZ2(rpad,tanth );
1106 clRphi=kACrphi-kBCrphi*rpad*tanth+kCCrphi*ratio*ratio;
1107 clZ=kACz-kBCz*rpad*tanth+kCCz*tanth*tanth;
1109 // temporary protection
1111 if(sigmaRphi < 0.) sigmaRphi=0.4e-3;
1112 if(sigmaZ < 0.) sigmaZ=0.4e-3;
1113 if(clRphi < 0.) clRphi=2.5e-3;
1114 if(clZ < 0.) clZ=2.5e-5;
1119 // smearing --> rotate to the 1 (13) or to the 25 (49) sector,
1120 // then the inaccuracy in a X-Y plane is only along Y (pad row)!
1122 Float_t xprim= tpcHit->X()*cph + tpcHit->Y()*sph;
1123 Float_t yprim=-tpcHit->X()*sph + tpcHit->Y()*cph;
1124 xyz[0]=gRandom->Gaus(yprim,TMath::Sqrt(sigmaRphi)); // y
1125 Float_t alpha=(isec < fTPCParam->GetNInnerSector()) ?
1126 fTPCParam->GetInnerAngle() : fTPCParam->GetOuterAngle();
1127 Float_t ymax=xprim*TMath::Tan(0.5*alpha);
1128 if (TMath::Abs(xyz[0])>ymax) xyz[0]=yprim;
1129 xyz[1]=gRandom->Gaus(tpcHit->Z(),TMath::Sqrt(sigmaZ)); // z
1130 if (TMath::Abs(xyz[1])>fTPCParam->GetZLength()) xyz[1]=tpcHit->Z();
1131 xyz[2]=sigmaRphi; // fSigmaY2
1132 xyz[3]=sigmaZ; // fSigmaZ2
1133 xyz[4]=tpcHit->fQ; // q
1135 AliTPCClustersRow *clrow=carray.GetRow(sector,tpcHit->fPadRow);
1136 if (!clrow) clrow=carray.CreateRow(sector,tpcHit->fPadRow);
1138 Int_t tracks[3]={tpcHit->Track(), -1, -1};
1139 AliTPCcluster cluster(tracks,xyz);
1141 clrow->InsertCluster(&cluster); nclusters++;
1143 tpcHit = (AliTPChit*)NextHit();
1146 } // end of loop over hits
1148 } // end of loop over tracks
1150 Int_t nrows=fTPCParam->GetNRow(isec);
1151 for (Int_t irow=0; irow<nrows; irow++) {
1152 AliTPCClustersRow *clrow=carray.GetRow(isec,irow);
1153 if (!clrow) continue;
1154 carray.StoreRow(isec,irow);
1155 carray.ClearRow(isec,irow);
1158 } // end of loop over sectors
1160 cerr<<"Number of made clusters : "<<nclusters<<" \n";
1161 carray.GetTree()->SetName(cname);
1162 carray.GetTree()->Write();
1164 savedir->cd(); //switch back to the input file
1166 } // end of function
1168 //_________________________________________________________________
1169 void AliTPC::Hits2ExactClustersSector(Int_t isec)
1171 //--------------------------------------------------------
1172 //calculate exact cross point of track and given pad row
1173 //resulting values are expressed in "digit" coordinata
1174 //--------------------------------------------------------
1176 //-----------------------------------------------------------------
1177 // Origin: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1178 //-----------------------------------------------------------------
1180 if (fClustersArray==0){
1184 TParticle *particle; // pointer to a given particle
1185 AliTPChit *tpcHit; // pointer to a sigle TPC hit
1186 // Int_t sector,nhits;
1188 const Int_t kcmaxhits=30000;
1189 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
1190 AliTPCFastVector & xxx = *xxxx;
1191 Int_t maxhits = kcmaxhits;
1192 //construct array for each padrow
1193 for (Int_t i=0; i<fTPCParam->GetNRow(isec);i++)
1194 fClustersArray->CreateRow(isec,i);
1196 //---------------------------------------------------------------
1197 // Get the access to the tracks
1198 //---------------------------------------------------------------
1200 TTree *tH = gAlice->TreeH();
1201 Stat_t ntracks = tH->GetEntries();
1202 Int_t npart = gAlice->GetNtrack();
1205 if (fHitType>1) branch = tH->GetBranch("TPC2");
1206 else branch = tH->GetBranch("TPC");
1208 //------------------------------------------------------------
1210 //------------------------------------------------------------
1212 for(Int_t track=0;track<ntracks;track++){
1213 Bool_t isInSector=kTRUE;
1215 isInSector = TrackInVolume(isec,track);
1216 if (!isInSector) continue;
1218 branch->GetEntry(track); // get next track
1220 // Get number of the TPC hits and a pointer
1224 Int_t currentIndex=0;
1225 Int_t lastrow=-1; //last writen row
1229 tpcHit = (AliTPChit*)FirstHit(-1);
1232 Int_t sector=tpcHit->fSector; // sector number
1234 tpcHit = (AliTPChit*) NextHit();
1238 ipart=tpcHit->Track();
1239 if (ipart<npart) particle=gAlice->Particle(ipart);
1243 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
1244 Int_t index[3]={1,isec,0};
1245 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
1246 if (currentrow<0) {tpcHit = (AliTPChit*)NextHit(); continue;}
1247 if (lastrow<0) lastrow=currentrow;
1248 if (currentrow==lastrow){
1249 if ( currentIndex>=maxhits){
1251 xxx.ResizeTo(4*maxhits);
1253 xxx(currentIndex*4)=x[0];
1254 xxx(currentIndex*4+1)=x[1];
1255 xxx(currentIndex*4+2)=x[2];
1256 xxx(currentIndex*4+3)=tpcHit->fQ;
1260 if (currentIndex>2){
1272 for (Int_t index=0;index<currentIndex;index++){
1274 x=x2=x3=x4=xxx(index*4);
1282 sumy+=xxx(index*4+1);
1283 sumxy+=xxx(index*4+1)*x;
1284 sumx2y+=xxx(index*4+1)*x2;
1285 sumz+=xxx(index*4+2);
1286 sumxz+=xxx(index*4+2)*x;
1287 sumx2z+=xxx(index*4+2)*x2;
1288 sumq+=xxx(index*4+3);
1290 Float_t centralPad = (fTPCParam->GetNPads(isec,lastrow)-1)/2;
1291 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
1292 sumx2*(sumx*sumx3-sumx2*sumx2);
1294 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
1295 sumx2*(sumxy*sumx3-sumx2y*sumx2);
1296 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
1297 sumx2*(sumxz*sumx3-sumx2z*sumx2);
1299 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
1300 sumx2*(sumx*sumx2y-sumx2*sumxy);
1301 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
1302 sumx2*(sumx*sumx2z-sumx2*sumxz);
1304 Float_t y=detay/det+centralPad;
1305 Float_t z=detaz/det;
1306 Float_t by=detby/det; //y angle
1307 Float_t bz=detbz/det; //z angle
1308 sumy/=Float_t(currentIndex);
1309 sumz/=Float_t(currentIndex);
1311 AliTPCClustersRow * row = (fClustersArray->GetRow(isec,lastrow));
1313 AliComplexCluster* cl = new((AliComplexCluster*)row->Append()) AliComplexCluster ;
1314 // AliComplexCluster cl;
1320 cl->fTracks[0]=ipart;
1324 } //end of calculating cluster for given row
1327 tpcHit = (AliTPChit*)NextHit();
1328 } // end of loop over hits
1329 } // end of loop over tracks
1330 //write padrows to tree
1331 for (Int_t ii=0; ii<fTPCParam->GetNRow(isec);ii++) {
1332 fClustersArray->StoreRow(isec,ii);
1333 fClustersArray->ClearRow(isec,ii);
1342 void AliTPC::SDigits2Digits2(Int_t eventnumber)
1344 //create digits from summable digits
1345 GenerNoise(500000); //create teble with noise
1348 sprintf(sname,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1349 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1351 //conect tree with sSDigits
1352 TTree *t = (TTree *)gDirectory->Get(sname);
1353 AliSimDigits digarr, *dummy=&digarr;
1354 t->GetBranch("Segment")->SetAddress(&dummy);
1355 Stat_t nentries = t->GetEntries();
1357 // set zero suppression
1359 fTPCParam->SetZeroSup(2);
1361 // get zero suppression
1363 Int_t zerosup = fTPCParam->GetZeroSup();
1365 //make tree with digits
1367 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1368 arr->SetClass("AliSimDigits");
1369 arr->Setup(fTPCParam);
1370 arr->MakeTree(fDigitsFile);
1372 AliTPCParam * par =fTPCParam;
1374 //Loop over segments of the TPC
1376 for (Int_t n=0; n<nentries; n++) {
1379 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
1380 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
1383 if (!IsSectorActive(sec)) continue;
1384 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1385 Int_t nrows = digrow->GetNRows();
1386 Int_t ncols = digrow->GetNCols();
1388 digrow->ExpandBuffer();
1389 digarr.ExpandBuffer();
1390 digrow->ExpandTrackBuffer();
1391 digarr.ExpandTrackBuffer();
1394 Short_t * pamp0 = digarr.GetDigits();
1395 Int_t * ptracks0 = digarr.GetTracks();
1396 Short_t * pamp1 = digrow->GetDigits();
1397 Int_t * ptracks1 = digrow->GetTracks();
1398 Int_t nelems =nrows*ncols;
1399 Int_t saturation = fTPCParam->GetADCSat();
1400 //use internal structure of the AliDigits - for speed reason
1401 //if you cahnge implementation
1402 //of the Alidigits - it must be rewriten -
1403 for (Int_t i= 0; i<nelems; i++){
1404 // Float_t q = *pamp0;
1405 //q/=16.; //conversion faktor
1406 //Float_t noise= GetNoise();
1408 //q= TMath::Nint(q);
1409 Float_t q = TMath::Nint(Float_t(*pamp0)/16.+GetNoise());
1411 if (q>saturation) q=saturation;
1413 //if (ptracks0[0]==0)
1416 ptracks1[0]=ptracks0[0];
1417 ptracks1[nelems]=ptracks0[nelems];
1418 ptracks1[2*nelems]=ptracks0[2*nelems];
1426 arr->StoreRow(sec,row);
1427 arr->ClearRow(sec,row);
1428 // cerr<<sec<<"\t"<<row<<"\n";
1435 arr->GetTree()->SetName(dname);
1436 arr->GetTree()->Write();
1439 //_________________________________________
1440 void AliTPC::Merge(TTree * intree, Int_t *mask, Int_t nin, Int_t outid)
1443 //intree - pointer to array of trees with s digits
1444 //mask - mask for each
1445 //nin - number of inputs
1446 //outid - event number of the output event
1448 //create digits from summable digits
1449 //conect tree with sSDigits
1452 AliSimDigits ** digarr = new AliSimDigits*[nin];
1453 for (Int_t i1=0;i1<nin; i1++){
1455 intree[i1].GetBranch("Segment")->SetAddress(&digarr[i1]);
1457 Stat_t nentries = intree[0].GetEntries();
1459 //make tree with digits
1461 sprintf(dname,"TreeD_%s_%d",fTPCParam->GetTitle(),outid);
1462 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1463 arr->SetClass("AliSimDigits");
1464 arr->Setup(fTPCParam);
1465 arr->MakeTree(fDigitsFile);
1467 // set zero suppression
1469 fTPCParam->SetZeroSup(2);
1471 // get zero suppression
1473 Int_t zerosup = fTPCParam->GetZeroSup();
1476 AliTPCParam * par =fTPCParam;
1478 //Loop over segments of the TPC
1479 for (Int_t n=0; n<nentries; n++) {
1481 for (Int_t i=0;i<nin; i++){
1482 //connect proper digits
1483 intree[i].GetEvent(n);
1484 digarr[i]->ExpandBuffer();
1485 digarr[i]->ExpandTrackBuffer();
1488 if (!par->AdjustSectorRow(digarr[0]->GetID(),sec,row)) {
1489 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr[0]->GetID()<<endl;
1493 AliSimDigits * digrow =(AliSimDigits*) arr->CreateRow(sec,row);
1494 Int_t nrows = digrow->GetNRows();
1495 Int_t ncols = digrow->GetNCols();
1497 digrow->ExpandBuffer();
1498 digrow->ExpandTrackBuffer();
1500 for (Int_t rows=0;rows<nrows; rows++){
1501 for (Int_t col=0;col<ncols; col++){
1503 Int_t label[1000]; // i hope no more than 300 events merged
1505 // looop over digits
1506 for (Int_t i=0;i<nin; i++){
1507 q += digarr[i]->GetDigitFast(rows,col);
1508 for (Int_t tr=0;tr<3;tr++) {
1509 Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
1511 label[labptr]=lab+mask[i]-2;
1517 q = gRandom->Gaus(q,fTPCParam->GetNoise()*fTPCParam->GetNoiseNormFac()*16);
1519 q/=16; //conversion faktor
1524 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat();
1525 digrow->SetDigitFast((Short_t)q,rows,col);
1526 for (Int_t tr=0;tr<3;tr++){
1528 ((AliSimDigits*)digrow)->SetTrackIDFast(label[tr],rows,col,tr);
1530 ((AliSimDigits*)digrow)->SetTrackIDFast(-1,rows,col,tr);
1535 arr->StoreRow(sec,row);
1537 arr->ClearRow(sec,row);
1542 arr->GetTree()->SetName(dname);
1543 arr->GetTree()->Write();
1549 /*_________________________________________
1550 void AliTPC::SDigits2Digits(Int_t eventnumber)
1554 cerr<<"Digitizing TPC...\n";
1556 Hits2Digits(eventnumber);
1561 // char treeName[100];
1563 // sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1565 // GetDigitsArray()->GetTree()->Write(treeName);
1568 //__________________________________________________________________
1569 void AliTPC::SetDefaults(){
1572 cerr<<"Setting default parameters...\n";
1574 // Set response functions
1576 AliTPCParamSR *param=(AliTPCParamSR*)gDirectory->Get("75x40_100x60");
1577 AliTPCPRF2D * prfinner = new AliTPCPRF2D;
1578 AliTPCPRF2D * prfouter1 = new AliTPCPRF2D;
1579 AliTPCPRF2D * prfouter2 = new AliTPCPRF2D;
1580 AliTPCRF1D * rf = new AliTPCRF1D(kTRUE);
1581 rf->SetGauss(param->GetZSigma(),param->GetZWidth(),1.);
1582 rf->SetOffset(3*param->GetZSigma());
1585 TDirectory *savedir=gDirectory;
1586 TFile *f=TFile::Open("$ALICE_ROOT/TPC/AliTPCprf2d.root");
1588 cerr<<"Can't open $ALICE_ROOT/TPC/AliTPCprf2d.root !\n" ;
1591 prfinner->Read("prf_07504_Gati_056068_d02");
1592 prfouter1->Read("prf_10006_Gati_047051_d03");
1593 prfouter2->Read("prf_15006_Gati_047051_d03");
1597 param->SetInnerPRF(prfinner);
1598 param->SetOuter1PRF(prfouter1);
1599 param->SetOuter2PRF(prfouter2);
1600 param->SetTimeRF(rf);
1610 //__________________________________________________________________
1611 void AliTPC::Hits2Digits(Int_t eventnumber)
1613 //----------------------------------------------------
1614 // Loop over all sectors for a single event
1615 //----------------------------------------------------
1618 if(fDefaults == 0) SetDefaults(); // check if the parameters are set
1619 GenerNoise(500000); //create teble with noise
1621 //setup TPCDigitsArray
1623 if(GetDigitsArray()) delete GetDigitsArray();
1625 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1626 arr->SetClass("AliSimDigits");
1627 arr->Setup(fTPCParam);
1628 arr->MakeTree(fDigitsFile);
1629 SetDigitsArray(arr);
1631 fDigitsSwitch=0; // standard digits
1633 cerr<<"Digitizing TPC -- normal digits...\n";
1635 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1641 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1643 GetDigitsArray()->GetTree()->SetName(treeName);
1644 GetDigitsArray()->GetTree()->Write();
1651 //__________________________________________________________________
1652 void AliTPC::Hits2SDigits2(Int_t eventnumber)
1655 //-----------------------------------------------------------
1656 // summable digits - 16 bit "ADC", no noise, no saturation
1657 //-----------------------------------------------------------
1659 //----------------------------------------------------
1660 // Loop over all sectors for a single event
1661 //----------------------------------------------------
1664 if(fDefaults == 0) SetDefaults();
1665 GenerNoise(500000); //create table with noise
1666 //setup TPCDigitsArray
1668 if(GetDigitsArray()) delete GetDigitsArray();
1670 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1671 arr->SetClass("AliSimDigits");
1672 arr->Setup(fTPCParam);
1673 arr->MakeTree(fDigitsFile);
1674 SetDigitsArray(arr);
1676 cerr<<"Digitizing TPC -- summable digits...\n";
1678 fDigitsSwitch=1; // summable digits
1680 // set zero suppression to "0"
1682 fTPCParam->SetZeroSup(0);
1684 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1691 sprintf(treeName,"TreeS_%s_%d",fTPCParam->GetTitle(),eventnumber);
1693 GetDigitsArray()->GetTree()->SetName(treeName);
1694 GetDigitsArray()->GetTree()->Write();
1700 //__________________________________________________________________
1701 void AliTPC::Hits2SDigits()
1704 //-----------------------------------------------------------
1705 // summable digits - 16 bit "ADC", no noise, no saturation
1706 //-----------------------------------------------------------
1708 //----------------------------------------------------
1709 // Loop over all sectors for a single event
1710 //----------------------------------------------------
1711 //MI change - for pp run
1712 Int_t eventnumber = gAlice->GetEvNumber();
1714 if(fDefaults == 0) SetDefaults();
1715 GenerNoise(500000); //create table with noise
1717 //setup TPCDigitsArray
1719 if(GetDigitsArray()) delete GetDigitsArray();
1721 AliTPCDigitsArray *arr = new AliTPCDigitsArray;
1722 arr->SetClass("AliSimDigits");
1723 arr->Setup(fTPCParam);
1724 arr->MakeTree(fDigitsFile);
1725 SetDigitsArray(arr);
1727 cerr<<"Digitizing TPC -- summable digits...\n";
1729 // fDigitsSwitch=1; // summable digits -for the moment direct
1731 for(Int_t isec=0;isec<fTPCParam->GetNSector();isec++) if (IsSectorActive(isec)) Hits2DigitsSector(isec);
1737 sprintf(treeName,"TreeD_%s_%d",fTPCParam->GetTitle(),eventnumber);
1739 GetDigitsArray()->GetTree()->SetName(treeName);
1740 GetDigitsArray()->GetTree()->Write();
1745 //_____________________________________________________________________________
1746 void AliTPC::Hits2DigitsSector(Int_t isec)
1748 //-------------------------------------------------------------------
1749 // TPC conversion from hits to digits.
1750 //-------------------------------------------------------------------
1752 //-----------------------------------------------------------------
1753 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1754 //-----------------------------------------------------------------
1756 //-------------------------------------------------------
1757 // Get the access to the track hits
1758 //-------------------------------------------------------
1760 // check if the parameters are set - important if one calls this method
1761 // directly, not from the Hits2Digits
1763 if(fDefaults == 0) SetDefaults();
1765 TTree *tH = gAlice->TreeH(); // pointer to the hits tree
1766 Stat_t ntracks = tH->GetEntries();
1770 //-------------------------------------------
1771 // Only if there are any tracks...
1772 //-------------------------------------------
1776 //printf("*** Processing sector number %d ***\n",isec);
1778 Int_t nrows =fTPCParam->GetNRow(isec);
1780 row= new TObjArray* [nrows];
1782 MakeSector(isec,nrows,tH,ntracks,row);
1784 //--------------------------------------------------------
1785 // Digitize this sector, row by row
1786 // row[i] is the pointer to the TObjArray of AliTPCFastVectors,
1787 // each one containing electrons accepted on this
1788 // row, assigned into tracks
1789 //--------------------------------------------------------
1793 if (fDigitsArray->GetTree()==0) fDigitsArray->MakeTree(fDigitsFile);
1795 for (i=0;i<nrows;i++){
1797 AliDigits * dig = fDigitsArray->CreateRow(isec,i);
1799 DigitizeRow(i,isec,row);
1801 fDigitsArray->StoreRow(isec,i);
1803 Int_t ndig = dig->GetDigitSize();
1804 if (gDebug > 10) printf("*** Sector, row, compressed digits %d %d %d ***\n",isec,i,ndig);
1806 fDigitsArray->ClearRow(isec,i);
1809 } // end of the sector digitization
1811 for(i=0;i<nrows;i++){
1816 delete [] row; // delete the array of pointers to TObjArray-s
1820 } // end of Hits2DigitsSector
1823 //_____________________________________________________________________________
1824 void AliTPC::DigitizeRow(Int_t irow,Int_t isec,TObjArray **rows)
1826 //-----------------------------------------------------------
1827 // Single row digitization, coupling from the neighbouring
1828 // rows taken into account
1829 //-----------------------------------------------------------
1831 //-----------------------------------------------------------------
1832 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1833 // Modified: Marian Ivanov GSI Darmstadt, m.ivanov@gsi.de
1834 //-----------------------------------------------------------------
1837 Float_t zerosup = fTPCParam->GetZeroSup();
1838 Int_t nrows =fTPCParam->GetNRow(isec);
1839 fCurrentIndex[1]= isec;
1842 Int_t nofPads = fTPCParam->GetNPads(isec,irow);
1843 Int_t nofTbins = fTPCParam->GetMaxTBin();
1844 Int_t indexRange[4];
1846 // Integrated signal for this row
1847 // and a single track signal
1850 AliTPCFastMatrix *m1 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // integrated
1851 AliTPCFastMatrix *m2 = new AliTPCFastMatrix(0,nofPads,0,nofTbins); // single
1853 AliTPCFastMatrix &total = *m1;
1855 // Array of pointers to the label-signal list
1857 Int_t nofDigits = nofPads*nofTbins; // number of digits for this row
1858 Float_t **pList = new Float_t* [nofDigits];
1862 for(lp=0;lp<nofDigits;lp++)pList[lp]=0; // set all pointers to NULL
1866 Int_t row1 = TMath::Max(irow-fTPCParam->GetNCrossRows(),0);
1867 Int_t row2 = TMath::Min(irow+fTPCParam->GetNCrossRows(),nrows-1);
1868 for (Int_t row= row1;row<=row2;row++){
1869 Int_t nTracks= rows[row]->GetEntries();
1870 for (i1=0;i1<nTracks;i1++){
1871 fCurrentIndex[2]= row;
1872 fCurrentIndex[3]=irow;
1874 m2->Zero(); // clear single track signal matrix
1875 Float_t trackLabel = GetSignal(rows[row],i1,m2,m1,indexRange);
1876 GetList(trackLabel,nofPads,m2,indexRange,pList);
1878 else GetSignal(rows[row],i1,0,m1,indexRange);
1884 AliDigits *dig = fDigitsArray->GetRow(isec,irow);
1886 Float_t fzerosup = zerosup+0.5;
1887 for(Int_t it=0;it<nofTbins;it++){
1888 Float_t *pq = &(total.UncheckedAt(0,it));
1889 for(Int_t ip=0;ip<nofPads;ip++){
1893 if(fDigitsSwitch == 0){
1895 if(q <=fzerosup) continue; // do not fill zeros
1897 if(q > fTPCParam->GetADCSat()) q = fTPCParam->GetADCSat(); // saturation
1902 if(q <= 0.) continue; // do not fill zeros
1903 if(q>2000.) q=2000.;
1909 // "real" signal or electronic noise (list = -1)?
1912 for(Int_t j1=0;j1<3;j1++){
1913 tracks[j1] = (pList[gi]) ?(Int_t)(*(pList[gi]+j1)) : -2;
1918 <A NAME="AliDigits"></A>
1919 using of AliDigits object
1922 dig->SetDigitFast((Short_t)q,it,ip);
1923 if (fDigitsArray->IsSimulated())
1925 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[0],it,ip,0);
1926 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[1],it,ip,1);
1927 ((AliSimDigits*)dig)->SetTrackIDFast(tracks[2],it,ip,2);
1931 } // end of loop over time buckets
1932 } // end of lop over pads
1935 // This row has been digitized, delete nonused stuff
1938 for(lp=0;lp<nofDigits;lp++){
1939 if(pList[lp]) delete [] pList[lp];
1948 } // end of DigitizeRow
1950 //_____________________________________________________________________________
1952 Float_t AliTPC::GetSignal(TObjArray *p1, Int_t ntr,
1953 AliTPCFastMatrix *m1, AliTPCFastMatrix *m2,Int_t *indexRange)
1956 //---------------------------------------------------------------
1957 // Calculates 2-D signal (pad,time) for a single track,
1958 // returns a pointer to the signal matrix and the track label
1959 // No digitization is performed at this level!!!
1960 //---------------------------------------------------------------
1962 //-----------------------------------------------------------------
1963 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
1964 // Modified: Marian Ivanov
1965 //-----------------------------------------------------------------
1967 AliTPCFastVector *tv;
1969 tv = (AliTPCFastVector*)p1->At(ntr); // pointer to a track
1970 AliTPCFastVector &v = *tv;
1972 Float_t label = v(0);
1973 Int_t centralPad = (fTPCParam->GetNPads(fCurrentIndex[1],fCurrentIndex[3])-1)/2;
1975 Int_t nElectrons = (tv->GetNrows()-1)/4;
1976 indexRange[0]=9999; // min pad
1977 indexRange[1]=-1; // max pad
1978 indexRange[2]=9999; //min time
1979 indexRange[3]=-1; // max time
1981 AliTPCFastMatrix &signal = *m1;
1982 AliTPCFastMatrix &total = *m2;
1984 // Loop over all electrons
1986 for(Int_t nel=0; nel<nElectrons; nel++){
1988 Float_t aval = v(idx+4);
1989 Float_t eltoadcfac=aval*fTPCParam->GetTotalNormFac();
1990 Float_t xyz[3]={v(idx+1),v(idx+2),v(idx+3)};
1991 Int_t n = ((AliTPCParamSR*)fTPCParam)->CalcResponseFast(xyz,fCurrentIndex,fCurrentIndex[3]);
1993 Int_t *index = fTPCParam->GetResBin(0);
1994 Float_t *weight = & (fTPCParam->GetResWeight(0));
1996 if (n>0) for (Int_t i =0; i<n; i++){
1997 Int_t pad=index[1]+centralPad; //in digit coordinates central pad has coordinate 0
2000 Int_t time=index[2];
2001 Float_t qweight = *(weight)*eltoadcfac;
2003 if (m1!=0) signal.UncheckedAt(pad,time)+=qweight;
2004 total.UncheckedAt(pad,time)+=qweight;
2005 if (indexRange[0]>pad) indexRange[0]=pad;
2006 if (indexRange[1]<pad) indexRange[1]=pad;
2007 if (indexRange[2]>time) indexRange[2]=time;
2008 if (indexRange[3]<time) indexRange[3]=time;
2015 } // end of loop over electrons
2017 return label; // returns track label when finished
2020 //_____________________________________________________________________________
2021 void AliTPC::GetList(Float_t label,Int_t np,AliTPCFastMatrix *m,
2022 Int_t *indexRange, Float_t **pList)
2024 //----------------------------------------------------------------------
2025 // Updates the list of tracks contributing to digits for a given row
2026 //----------------------------------------------------------------------
2028 //-----------------------------------------------------------------
2029 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2030 //-----------------------------------------------------------------
2032 AliTPCFastMatrix &signal = *m;
2034 // lop over nonzero digits
2036 for(Int_t it=indexRange[2];it<indexRange[3]+1;it++){
2037 for(Int_t ip=indexRange[0];ip<indexRange[1]+1;ip++){
2040 // accept only the contribution larger than 500 electrons (1/2 s_noise)
2042 if(signal(ip,it)<0.5) continue;
2045 Int_t globalIndex = it*np+ip; // globalIndex starts from 0!
2047 if(!pList[globalIndex]){
2050 // Create new list (6 elements - 3 signals and 3 labels),
2053 pList[globalIndex] = new Float_t [6];
2057 *pList[globalIndex] = -1.;
2058 *(pList[globalIndex]+1) = -1.;
2059 *(pList[globalIndex]+2) = -1.;
2060 *(pList[globalIndex]+3) = -1.;
2061 *(pList[globalIndex]+4) = -1.;
2062 *(pList[globalIndex]+5) = -1.;
2065 *pList[globalIndex] = label;
2066 *(pList[globalIndex]+3) = signal(ip,it);
2070 // check the signal magnitude
2072 Float_t highest = *(pList[globalIndex]+3);
2073 Float_t middle = *(pList[globalIndex]+4);
2074 Float_t lowest = *(pList[globalIndex]+5);
2077 // compare the new signal with already existing list
2080 if(signal(ip,it)<lowest) continue; // neglect this track
2084 if (signal(ip,it)>highest){
2085 *(pList[globalIndex]+5) = middle;
2086 *(pList[globalIndex]+4) = highest;
2087 *(pList[globalIndex]+3) = signal(ip,it);
2089 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2090 *(pList[globalIndex]+1) = *pList[globalIndex];
2091 *pList[globalIndex] = label;
2093 else if (signal(ip,it)>middle){
2094 *(pList[globalIndex]+5) = middle;
2095 *(pList[globalIndex]+4) = signal(ip,it);
2097 *(pList[globalIndex]+2) = *(pList[globalIndex]+1);
2098 *(pList[globalIndex]+1) = label;
2101 *(pList[globalIndex]+5) = signal(ip,it);
2102 *(pList[globalIndex]+2) = label;
2106 } // end of loop over pads
2107 } // end of loop over time bins
2112 //___________________________________________________________________
2113 void AliTPC::MakeSector(Int_t isec,Int_t nrows,TTree *TH,
2114 Stat_t ntracks,TObjArray **row)
2117 //-----------------------------------------------------------------
2118 // Prepares the sector digitization, creates the vectors of
2119 // tracks for each row of this sector. The track vector
2120 // contains the track label and the position of electrons.
2121 //-----------------------------------------------------------------
2123 //-----------------------------------------------------------------
2124 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2125 //-----------------------------------------------------------------
2127 Float_t gasgain = fTPCParam->GetGasGain();
2131 AliTPChit *tpcHit; // pointer to a sigle TPC hit
2134 if (fHitType>1) branch = TH->GetBranch("TPC2");
2135 else branch = TH->GetBranch("TPC");
2138 //----------------------------------------------
2139 // Create TObjArray-s, one for each row,
2140 // each TObjArray will store the AliTPCFastVectors
2141 // of electrons, one AliTPCFastVectors per each track.
2142 //----------------------------------------------
2144 Int_t *nofElectrons = new Int_t [nrows]; // electron counter for each row
2145 AliTPCFastVector **tracks = new AliTPCFastVector* [nrows]; //pointers to the track vectors
2147 for(i=0; i<nrows; i++){
2148 row[i] = new TObjArray;
2155 //--------------------------------------------------------------------
2156 // Loop over tracks, the "track" contains the full history
2157 //--------------------------------------------------------------------
2159 Int_t previousTrack,currentTrack;
2160 previousTrack = -1; // nothing to store so far!
2162 for(Int_t track=0;track<ntracks;track++){
2163 Bool_t isInSector=kTRUE;
2165 isInSector = TrackInVolume(isec,track);
2166 if (!isInSector) continue;
2168 branch->GetEntry(track); // get next track
2172 tpcHit = (AliTPChit*)FirstHit(-1);
2174 //--------------------------------------------------------------
2176 //--------------------------------------------------------------
2181 Int_t sector=tpcHit->fSector; // sector number
2183 tpcHit = (AliTPChit*) NextHit();
2187 currentTrack = tpcHit->Track(); // track number
2190 if(currentTrack != previousTrack){
2192 // store already filled fTrack
2194 for(i=0;i<nrows;i++){
2195 if(previousTrack != -1){
2196 if(nofElectrons[i]>0){
2197 AliTPCFastVector &v = *tracks[i];
2198 v(0) = previousTrack;
2199 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2200 row[i]->Add(tracks[i]);
2203 delete tracks[i]; // delete empty AliTPCFastVector
2209 tracks[i] = new AliTPCFastVector(481); // AliTPCFastVectors for the next fTrack
2211 } // end of loop over rows
2213 previousTrack=currentTrack; // update track label
2216 Int_t qI = (Int_t) (tpcHit->fQ); // energy loss (number of electrons)
2218 //---------------------------------------------------
2219 // Calculate the electron attachment probability
2220 //---------------------------------------------------
2223 Float_t time = 1.e6*(fTPCParam->GetZLength()-TMath::Abs(tpcHit->Z()))
2224 /fTPCParam->GetDriftV();
2226 Float_t attProb = fTPCParam->GetAttCoef()*
2227 fTPCParam->GetOxyCont()*time; // fraction!
2229 //-----------------------------------------------
2230 // Loop over electrons
2231 //-----------------------------------------------
2234 for(Int_t nel=0;nel<qI;nel++){
2235 // skip if electron lost due to the attachment
2236 if((gRandom->Rndm(0)) < attProb) continue; // electron lost!
2241 // protection for the nonphysical avalanche size (10**6 maximum)
2243 Double_t rn=TMath::Max(gRandom->Rndm(0),1.93e-22);
2244 xyz[3]= (Float_t) (-gasgain*TMath::Log(rn));
2247 TransportElectron(xyz,index); //MI change -august
2249 fTPCParam->GetPadRow(xyz,index); //MI change august
2250 rowNumber = index[2];
2251 //transform position to local digit coordinates
2252 //relative to nearest pad row
2253 if ((rowNumber<0)||rowNumber>=fTPCParam->GetNRow(isec)) continue;
2255 if (isec <fTPCParam->GetNInnerSector()) {
2256 x1 = xyz[1]*fTPCParam->GetInnerPadPitchWidth();
2257 y1 = fTPCParam->GetYInner(rowNumber);
2260 x1=xyz[1]*fTPCParam->GetOuterPadPitchWidth();
2261 y1 = fTPCParam->GetYOuter(rowNumber);
2271 nofElectrons[rowNumber]++;
2272 //----------------------------------
2273 // Expand vector if necessary
2274 //----------------------------------
2275 if(nofElectrons[rowNumber]>120){
2276 Int_t range = tracks[rowNumber]->GetNrows();
2277 if((nofElectrons[rowNumber])>(range-1)/4){
2279 tracks[rowNumber]->ResizeTo(range+400); // Add 100 electrons
2283 AliTPCFastVector &v = *tracks[rowNumber];
2284 Int_t idx = 4*nofElectrons[rowNumber]-3;
2285 Real_t * position = &(((AliTPCFastVector&)v).UncheckedAt(idx)); //make code faster
2286 memcpy(position,xyz,4*sizeof(Float_t));
2288 } // end of loop over electrons
2290 tpcHit = (AliTPChit*)NextHit();
2292 } // end of loop over hits
2293 } // end of loop over tracks
2296 // store remaining track (the last one) if not empty
2299 for(i=0;i<nrows;i++){
2300 if(nofElectrons[i]>0){
2301 AliTPCFastVector &v = *tracks[i];
2302 v(0) = previousTrack;
2303 tracks[i]->ResizeTo(4*nofElectrons[i]+1); // shrink if necessary
2304 row[i]->Add(tracks[i]);
2313 delete [] nofElectrons;
2316 } // end of MakeSector
2319 //_____________________________________________________________________________
2323 // Initialise TPC detector after definition of geometry
2328 printf("\n%s: ",ClassName());
2329 for(i=0;i<35;i++) printf("*");
2330 printf(" TPC_INIT ");
2331 for(i=0;i<35;i++) printf("*");
2332 printf("\n%s: ",ClassName());
2334 for(i=0;i<80;i++) printf("*");
2339 //_____________________________________________________________________________
2340 void AliTPC::MakeBranch(Option_t* option, const char *file)
2343 // Create Tree branches for the TPC.
2345 Int_t buffersize = 4000;
2346 char branchname[10];
2347 sprintf(branchname,"%s",GetName());
2349 AliDetector::MakeBranch(option,file);
2351 const char *d = strstr(option,"D");
2353 if (fDigits && gAlice->TreeD() && d) {
2354 MakeBranchInTree(gAlice->TreeD(),
2355 branchname, &fDigits, buffersize, file);
2358 if (fHitType>1) MakeBranch2(option,file); // MI change 14.09.2000
2361 //_____________________________________________________________________________
2362 void AliTPC::ResetDigits()
2365 // Reset number of digits and the digits array for this detector
2368 if (fDigits) fDigits->Clear();
2371 //_____________________________________________________________________________
2372 void AliTPC::SetSecAL(Int_t sec)
2374 //---------------------------------------------------
2375 // Activate/deactivate selection for lower sectors
2376 //---------------------------------------------------
2378 //-----------------------------------------------------------------
2379 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2380 //-----------------------------------------------------------------
2385 //_____________________________________________________________________________
2386 void AliTPC::SetSecAU(Int_t sec)
2388 //----------------------------------------------------
2389 // Activate/deactivate selection for upper sectors
2390 //---------------------------------------------------
2392 //-----------------------------------------------------------------
2393 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2394 //-----------------------------------------------------------------
2399 //_____________________________________________________________________________
2400 void AliTPC::SetSecLows(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6)
2402 //----------------------------------------
2403 // Select active lower sectors
2404 //----------------------------------------
2406 //-----------------------------------------------------------------
2407 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2408 //-----------------------------------------------------------------
2418 //_____________________________________________________________________________
2419 void AliTPC::SetSecUps(Int_t s1,Int_t s2,Int_t s3,Int_t s4,Int_t s5, Int_t s6,
2420 Int_t s7, Int_t s8 ,Int_t s9 ,Int_t s10,
2421 Int_t s11 , Int_t s12)
2423 //--------------------------------
2424 // Select active upper sectors
2425 //--------------------------------
2427 //-----------------------------------------------------------------
2428 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2429 //-----------------------------------------------------------------
2445 //_____________________________________________________________________________
2446 void AliTPC::SetSens(Int_t sens)
2449 //-------------------------------------------------------------
2450 // Activates/deactivates the sensitive strips at the center of
2451 // the pad row -- this is for the space-point resolution calculations
2452 //-------------------------------------------------------------
2454 //-----------------------------------------------------------------
2455 // Origin: Marek Kowalski IFJ, Krakow, Marek.Kowalski@ifj.edu.pl
2456 //-----------------------------------------------------------------
2462 void AliTPC::SetSide(Float_t side=0.)
2464 // choice of the TPC side
2469 //____________________________________________________________________________
2470 void AliTPC::SetGasMixt(Int_t nc,Int_t c1,Int_t c2,Int_t c3,Float_t p1,
2471 Float_t p2,Float_t p3)
2474 // gax mixture definition
2488 //_____________________________________________________________________________
2490 void AliTPC::TransportElectron(Float_t *xyz, Int_t *index)
2493 // electron transport taking into account:
2495 // 2.ExB at the wires
2496 // 3. nonisochronity
2498 // xyz and index must be already transformed to system 1
2501 fTPCParam->Transform1to2(xyz,index);
2504 Float_t driftl=xyz[2];
2505 if(driftl<0.01) driftl=0.01;
2506 driftl=TMath::Sqrt(driftl);
2507 Float_t sigT = driftl*(fTPCParam->GetDiffT());
2508 Float_t sigL = driftl*(fTPCParam->GetDiffL());
2509 xyz[0]=gRandom->Gaus(xyz[0],sigT);
2510 xyz[1]=gRandom->Gaus(xyz[1],sigT);
2511 xyz[2]=gRandom->Gaus(xyz[2],sigL);
2515 if (fTPCParam->GetMWPCReadout()==kTRUE){
2517 fTPCParam->Transform2to2NearestWire(xyz,index);
2518 Float_t dx=xyz[0]-x1;
2519 xyz[1]+=dx*(fTPCParam->GetOmegaTau());
2521 //add nonisochronity (not implemented yet)
2525 ClassImp(AliTPCdigit)
2527 //_____________________________________________________________________________
2528 AliTPCdigit::AliTPCdigit(Int_t *tracks, Int_t *digits):
2532 // Creates a TPC digit object
2534 fSector = digits[0];
2535 fPadRow = digits[1];
2538 fSignal = digits[4];
2544 //_____________________________________________________________________________
2545 AliTPChit::AliTPChit(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
2549 // Creates a TPC hit object
2560 //________________________________________________________________________
2561 // Additional code because of the AliTPCTrackHitsV2
2563 void AliTPC::MakeBranch2(Option_t *option,const char *file)
2566 // Create a new branch in the current Root Tree
2567 // The branch of fHits is automatically split
2568 // MI change 14.09.2000
2569 if (fHitType<2) return;
2570 char branchname[10];
2571 sprintf(branchname,"%s2",GetName());
2573 // Get the pointer to the header
2574 const char *cH = strstr(option,"H");
2576 if (fTrackHits && gAlice->TreeH() && cH && fHitType&4) {
2577 // AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2578 // gAlice->TreeH(),fBufferSize,99);
2579 //TBranch * branch =
2580 gAlice->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2583 // gAlice->TreeH()->GetListOfBranches()->Add(branch);
2585 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2586 const char folder [] = "RunMC/Event/Data";
2588 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2589 Publish(folder,&fTrackHits,branchname);
2591 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2592 TDirectory *wd = gDirectory;
2594 TIter next( b->GetListOfBranches());
2595 while ((b=(TBranch*)next())) {
2600 cout << "Diverting branch " << branchname << " to file " << file << endl;
2604 if (fTrackHitsOld && gAlice->TreeH() && cH && fHitType&2) {
2605 AliObjectBranch * branch = new AliObjectBranch(branchname,"AliTPCTrackHits",&fTrackHitsOld,
2606 gAlice->TreeH(),fBufferSize,99);
2607 //TBranch * branch = gAlice->TreeH()->Branch(branchname,"AliTPCTrackHitsV2",&fTrackHits,
2610 gAlice->TreeH()->GetListOfBranches()->Add(branch);
2612 printf("* AliDetector::MakeBranch * Making Branch %s for trackhits\n",branchname);
2613 const char folder [] = "RunMC/Event/Data";
2615 printf("%15s: Publishing %s to %s\n",ClassName(),branchname,folder);
2616 Publish(folder,&fTrackHitsOld,branchname);
2618 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
2619 TDirectory *wd = gDirectory;
2621 TIter next( b->GetListOfBranches());
2622 while ((b=(TBranch*)next())) {
2627 cout << "Diverting branch " << branchname << " to file " << file << endl;
2632 void AliTPC::SetTreeAddress()
2634 if (fHitType<=1) AliDetector::SetTreeAddress();
2635 if (fHitType>1) SetTreeAddress2();
2638 void AliTPC::SetTreeAddress2()
2641 // Set branch address for the TrackHits Tree
2644 char branchname[20];
2645 sprintf(branchname,"%s2",GetName());
2647 // Branch address for hit tree
2648 TTree *treeH = gAlice->TreeH();
2649 if ((treeH)&&(fHitType&4)) {
2650 branch = treeH->GetBranch(branchname);
2651 if (branch) branch->SetAddress(&fTrackHits);
2653 if ((treeH)&&(fHitType&2)) {
2654 branch = treeH->GetBranch(branchname);
2655 if (branch) branch->SetAddress(&fTrackHitsOld);
2660 void AliTPC::FinishPrimary()
2662 if (fTrackHits &&fHitType&4) fTrackHits->FlushHitStack();
2663 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->FlushHitStack();
2667 void AliTPC::AddHit2(Int_t track, Int_t *vol, Float_t *hits)
2670 // add hit to the list
2673 int primary = gAlice->GetPrimary(track);
2674 gAlice->Particle(primary)->SetBit(kKeepBit);
2678 gAlice->FlagTrack(track);
2680 //AliTPChit *hit = (AliTPChit*)fHits->UncheckedAt(fNhits-1);
2681 //if (hit->fTrack!=rtrack)
2682 // cout<<"bad track number\n";
2683 if (fTrackHits && fHitType&4)
2684 fTrackHits->AddHitKartez(vol[0],rtrack, hits[0],
2685 hits[1],hits[2],(Int_t)hits[3]);
2686 if (fTrackHitsOld &&fHitType&2 )
2687 fTrackHitsOld->AddHitKartez(vol[0],rtrack, hits[0],
2688 hits[1],hits[2],(Int_t)hits[3]);
2692 void AliTPC::ResetHits()
2694 if (fHitType&1) AliDetector::ResetHits();
2695 if (fHitType>1) ResetHits2();
2698 void AliTPC::ResetHits2()
2702 if (fTrackHits && fHitType&4) fTrackHits->Clear();
2703 if (fTrackHitsOld && fHitType&2) fTrackHitsOld->Clear();
2707 AliHit* AliTPC::FirstHit(Int_t track)
2709 if (fHitType>1) return FirstHit2(track);
2710 return AliDetector::FirstHit(track);
2712 AliHit* AliTPC::NextHit()
2714 if (fHitType>1) return NextHit2();
2716 return AliDetector::NextHit();
2719 AliHit* AliTPC::FirstHit2(Int_t track)
2722 // Initialise the hit iterator
2723 // Return the address of the first hit for track
2724 // If track>=0 the track is read from disk
2725 // while if track<0 the first hit of the current
2726 // track is returned
2729 gAlice->ResetHits();
2730 gAlice->TreeH()->GetEvent(track);
2733 if (fTrackHits && fHitType&4) {
2734 fTrackHits->First();
2735 return fTrackHits->GetHit();
2737 if (fTrackHitsOld && fHitType&2) {
2738 fTrackHitsOld->First();
2739 return fTrackHitsOld->GetHit();
2745 AliHit* AliTPC::NextHit2()
2748 //Return the next hit for the current track
2751 if (fTrackHitsOld && fHitType&2) {
2752 fTrackHitsOld->Next();
2753 return fTrackHitsOld->GetHit();
2757 return fTrackHits->GetHit();
2763 void AliTPC::LoadPoints(Int_t)
2767 /* if(fHitType==1) return AliDetector::LoadPoints(a);
2770 if(fHitType==1) AliDetector::LoadPoints(a);
2771 else LoadPoints2(a);
2778 void AliTPC::RemapTrackHitIDs(Int_t *map)
2780 if (!fTrackHits) return;
2782 if (fTrackHitsOld && fHitType&2){
2783 AliObjectArray * arr = fTrackHitsOld->fTrackHitsInfo;
2784 for (UInt_t i=0;i<arr->GetSize();i++){
2785 AliTrackHitsInfo * info = (AliTrackHitsInfo *)(arr->At(i));
2786 info->fTrackID = map[info->fTrackID];
2789 if (fTrackHitsOld && fHitType&4){
2790 TClonesArray * arr = fTrackHits->GetArray();;
2791 for (Int_t i=0;i<arr->GetEntriesFast();i++){
2792 AliTrackHitsParamV2 * info = (AliTrackHitsParamV2 *)(arr->At(i));
2793 info->fTrackID = map[info->fTrackID];
2798 Bool_t AliTPC::TrackInVolume(Int_t id,Int_t track)
2800 //return bool information - is track in given volume
2801 //load only part of the track information
2802 //return true if current track is in volume
2805 if (fTrackHitsOld && fHitType&2) {
2806 TBranch * br = gAlice->TreeH()->GetBranch("fTrackHitsInfo");
2807 br->GetEvent(track);
2808 AliObjectArray * ar = fTrackHitsOld->fTrackHitsInfo;
2809 for (UInt_t j=0;j<ar->GetSize();j++){
2810 if ( ((AliTrackHitsInfo*)ar->At(j))->fVolumeID==id) return kTRUE;
2814 if (fTrackHits && fHitType&4) {
2815 TBranch * br1 = gAlice->TreeH()->GetBranch("fVolumes");
2816 TBranch * br2 = gAlice->TreeH()->GetBranch("fNVolumes");
2817 br2->GetEvent(track);
2818 br1->GetEvent(track);
2819 Int_t *volumes = fTrackHits->GetVolumes();
2820 Int_t nvolumes = fTrackHits->GetNVolumes();
2821 if (!volumes && nvolumes>0) {
2822 printf("Problematic track\t%d\t%d",track,nvolumes);
2825 for (Int_t j=0;j<nvolumes; j++)
2826 if (volumes[j]==id) return kTRUE;
2830 TBranch * br = gAlice->TreeH()->GetBranch("fSector");
2831 br->GetEvent(track);
2832 for (Int_t j=0;j<fHits->GetEntriesFast();j++){
2833 if ( ((AliTPChit*)fHits->At(j))->fSector==id) return kTRUE;
2840 //_____________________________________________________________________________
2841 void AliTPC::LoadPoints2(Int_t)
2844 // Store x, y, z of all hits in memory
2846 if (fTrackHits == 0 && fTrackHitsOld==0) return;
2849 if (fHitType&4) nhits = fTrackHits->GetEntriesFast();
2850 if (fHitType&2) nhits = fTrackHitsOld->GetEntriesFast();
2852 if (nhits == 0) return;
2853 Int_t tracks = gAlice->GetNtrack();
2854 if (fPoints == 0) fPoints = new TObjArray(tracks);
2857 Int_t *ntrk=new Int_t[tracks];
2858 Int_t *limi=new Int_t[tracks];
2859 Float_t **coor=new Float_t*[tracks];
2860 for(Int_t i=0;i<tracks;i++) {
2866 AliPoints *points = 0;
2869 Int_t chunk=nhits/4+1;
2871 // Loop over all the hits and store their position
2873 ahit = FirstHit2(-1);
2875 trk=ahit->GetTrack();
2876 if(ntrk[trk]==limi[trk]) {
2878 // Initialise a new track
2879 fp=new Float_t[3*(limi[trk]+chunk)];
2881 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2882 delete [] coor[trk];
2889 fp[3*ntrk[trk] ] = ahit->X();
2890 fp[3*ntrk[trk]+1] = ahit->Y();
2891 fp[3*ntrk[trk]+2] = ahit->Z();
2899 for(trk=0; trk<tracks; ++trk) {
2901 points = new AliPoints();
2902 points->SetMarkerColor(GetMarkerColor());
2903 points->SetMarkerSize(GetMarkerSize());
2904 points->SetDetector(this);
2905 points->SetParticle(trk);
2906 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
2907 fPoints->AddAt(points,trk);
2908 delete [] coor[trk];
2918 //_____________________________________________________________________________
2919 void AliTPC::LoadPoints3(Int_t)
2922 // Store x, y, z of all hits in memory
2923 // - only intersection point with pad row
2924 if (fTrackHits == 0) return;
2926 Int_t nhits = fTrackHits->GetEntriesFast();
2927 if (nhits == 0) return;
2928 Int_t tracks = gAlice->GetNtrack();
2929 if (fPoints == 0) fPoints = new TObjArray(2*tracks);
2930 fPoints->Expand(2*tracks);
2933 Int_t *ntrk=new Int_t[tracks];
2934 Int_t *limi=new Int_t[tracks];
2935 Float_t **coor=new Float_t*[tracks];
2936 for(Int_t i=0;i<tracks;i++) {
2942 AliPoints *points = 0;
2945 Int_t chunk=nhits/4+1;
2947 // Loop over all the hits and store their position
2949 ahit = FirstHit2(-1);
2950 //for (Int_t hit=0;hit<nhits;hit++) {
2954 // ahit = (AliHit*)fHits->UncheckedAt(hit);
2955 trk=ahit->GetTrack();
2956 Float_t x[3]={ahit->X(),ahit->Y(),ahit->Z()};
2957 Int_t index[3]={1,((AliTPChit*)ahit)->fSector,0};
2958 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
2959 if (currentrow!=lastrow){
2960 lastrow = currentrow;
2961 //later calculate intersection point
2962 if(ntrk[trk]==limi[trk]) {
2964 // Initialise a new track
2965 fp=new Float_t[3*(limi[trk]+chunk)];
2967 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
2968 delete [] coor[trk];
2975 fp[3*ntrk[trk] ] = ahit->X();
2976 fp[3*ntrk[trk]+1] = ahit->Y();
2977 fp[3*ntrk[trk]+2] = ahit->Z();
2984 for(trk=0; trk<tracks; ++trk) {
2986 points = new AliPoints();
2987 points->SetMarkerColor(GetMarkerColor()+1);
2988 points->SetMarkerStyle(5);
2989 points->SetMarkerSize(0.2);
2990 points->SetDetector(this);
2991 points->SetParticle(trk);
2992 // points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle()20);
2993 points->SetPolyMarker(ntrk[trk],coor[trk],30);
2994 fPoints->AddAt(points,tracks+trk);
2995 delete [] coor[trk];
3006 void AliTPC::FindTrackHitsIntersection(TClonesArray * arr)
3010 //fill clones array with intersection of current point with the
3015 const Int_t kcmaxhits=30000;
3016 AliTPCFastVector * xxxx = new AliTPCFastVector(kcmaxhits*4);
3017 AliTPCFastVector & xxx = *xxxx;
3018 Int_t maxhits = kcmaxhits;
3021 AliTPChit * tpcHit=0;
3022 tpcHit = (AliTPChit*)FirstHit2(-1);
3023 Int_t currentIndex=0;
3024 Int_t lastrow=-1; //last writen row
3027 if (tpcHit==0) continue;
3028 sector=tpcHit->fSector; // sector number
3029 ipart=tpcHit->Track();
3033 Float_t x[3]={tpcHit->X(),tpcHit->Y(),tpcHit->Z()};
3034 Int_t index[3]={1,sector,0};
3035 Int_t currentrow = fTPCParam->GetPadRow(x,index) ;
3036 if (currentrow<0) continue;
3037 if (lastrow<0) lastrow=currentrow;
3038 if (currentrow==lastrow){
3039 if ( currentIndex>=maxhits){
3041 xxx.ResizeTo(4*maxhits);
3043 xxx(currentIndex*4)=x[0];
3044 xxx(currentIndex*4+1)=x[1];
3045 xxx(currentIndex*4+2)=x[2];
3046 xxx(currentIndex*4+3)=tpcHit->fQ;
3050 if (currentIndex>2){
3062 for (Int_t index=0;index<currentIndex;index++){
3064 x=x2=x3=x4=xxx(index*4);
3072 sumy+=xxx(index*4+1);
3073 sumxy+=xxx(index*4+1)*x;
3074 sumx2y+=xxx(index*4+1)*x2;
3075 sumz+=xxx(index*4+2);
3076 sumxz+=xxx(index*4+2)*x;
3077 sumx2z+=xxx(index*4+2)*x2;
3078 sumq+=xxx(index*4+3);
3080 Float_t centralPad = (fTPCParam->GetNPads(sector,lastrow)-1)/2;
3081 Float_t det=currentIndex*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumx*sumx4-sumx2*sumx3)+
3082 sumx2*(sumx*sumx3-sumx2*sumx2);
3084 Float_t detay=sumy*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxy*sumx4-sumx2y*sumx3)+
3085 sumx2*(sumxy*sumx3-sumx2y*sumx2);
3086 Float_t detaz=sumz*(sumx2*sumx4-sumx3*sumx3)-sumx*(sumxz*sumx4-sumx2z*sumx3)+
3087 sumx2*(sumxz*sumx3-sumx2z*sumx2);
3089 Float_t detby=currentIndex*(sumxy*sumx4-sumx2y*sumx3)-sumy*(sumx*sumx4-sumx2*sumx3)+
3090 sumx2*(sumx*sumx2y-sumx2*sumxy);
3091 Float_t detbz=currentIndex*(sumxz*sumx4-sumx2z*sumx3)-sumz*(sumx*sumx4-sumx2*sumx3)+
3092 sumx2*(sumx*sumx2z-sumx2*sumxz);
3094 Float_t y=detay/det+centralPad;
3095 Float_t z=detaz/det;
3096 Float_t by=detby/det; //y angle
3097 Float_t bz=detbz/det; //z angle
3098 sumy/=Float_t(currentIndex);
3099 sumz/=Float_t(currentIndex);
3101 AliComplexCluster cl;
3107 cl.fTracks[0]=ipart;
3109 AliTPCClustersRow * row = (fClustersArray->GetRow(sector,lastrow));
3110 if (row!=0) row->InsertCluster(&cl);
3113 } //end of calculating cluster for given row
3115 } // end of loop over hits
3119 //_______________________________________________________________________________
3120 void AliTPC::Digits2Reco(Int_t firstevent,Int_t lastevent)
3122 // produces rec points from digits and writes them on the root file
3123 // AliTPCclusters.root
3125 TDirectory *cwd = gDirectory;
3128 AliTPCParam *dig=(AliTPCParam *)gDirectory->Get("75x40_100x60");
3130 cout<<"AliTPC::Digits2Reco: TPC parameteres have been set"<<endl;
3132 if(!gSystem->Getenv("CONFIG_FILE")){
3133 out=TFile::Open("AliTPCclusters.root","recreate");
3139 tmp1=gSystem->Getenv("CONFIG_FILE_PREFIX");
3140 tmp2=gSystem->Getenv("CONFIG_OUTDIR");
3141 sprintf(tmp3,"%s%s/AliTPCclusters.root",tmp1,tmp2);
3142 out=TFile::Open(tmp3,"recreate");
3146 cout<<"AliTPC::Digits2Reco - determination of rec points begins"<<endl;
3149 for(Int_t iev=firstevent;iev<lastevent+1;iev++){
3151 printf("Processing event %d\n",iev);
3152 Digits2Clusters(out,iev);
3154 cout<<"AliTPC::Digits2Reco - determination of rec points ended"<<endl;